Actual source code: ex44.c
2: #include <petscvec.h>
4: static char help[] = "Tests vecScatter Sequential to Sequential for (CUDA) vectors\n\
5: -m # : the size of the vectors\n \
6: -n # : the number of indices (with n<=m)\n \
7: -toFirst # : the starting index of the output vector for strided scatters\n \
8: -toStep # : the step size into the output vector for strided scatters\n \
9: -fromFirst # : the starting index of the input vector for strided scatters\n\
10: -fromStep # : the step size into the input vector for strided scatters\n\n";
12: int main(int argc, char *argv[])
13: {
14: Vec X, Y;
15: PetscInt m, n, i, n1, n2;
16: PetscInt toFirst, toStep, fromFirst, fromStep;
17: PetscInt *idx, *idy;
18: PetscBool flg;
19: IS toISStrided, fromISStrided, toISGeneral, fromISGeneral;
20: VecScatter vscatSStoSS, vscatSStoSG, vscatSGtoSS, vscatSGtoSG;
21: ScatterMode mode;
22: InsertMode addv;
24: PetscFunctionBeginUser;
25: PetscCall(PetscInitialize(&argc, &argv, 0, help));
26: PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, &flg));
27: if (!flg) m = 100;
29: PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, &flg));
30: if (!flg) n = 30;
32: PetscCall(PetscOptionsGetInt(NULL, NULL, "-toFirst", &toFirst, &flg));
33: if (!flg) toFirst = 3;
35: PetscCall(PetscOptionsGetInt(NULL, NULL, "-toStep", &toStep, &flg));
36: if (!flg) toStep = 3;
38: PetscCall(PetscOptionsGetInt(NULL, NULL, "-fromFirst", &fromFirst, &flg));
39: if (!flg) fromFirst = 2;
41: PetscCall(PetscOptionsGetInt(NULL, NULL, "-fromStep", &fromStep, &flg));
42: if (!flg) fromStep = 2;
44: if (n > m) {
45: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "The vector sizes are %" PetscInt_FMT ". The number of elements being scattered is %" PetscInt_FMT "\n", m, n));
46: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Adjust the parameters such that m>=n\n"));
47: } else if (toFirst + (n - 1) * toStep >= m) {
48: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "The vector sizes are %" PetscInt_FMT ". The number of elements being scattered is %" PetscInt_FMT "\n", m, n));
49: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "For the Strided Scatter, toFirst=%" PetscInt_FMT " and toStep=%" PetscInt_FMT ".\n", toFirst, toStep));
50: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "This produces an index (toFirst+(n-1)*toStep)>=m\n"));
51: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Adjust the parameterrs accordingly with -m, -n, -toFirst, or -toStep\n"));
52: } else if (fromFirst + (n - 1) * fromStep >= m) {
53: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "The vector sizes are %" PetscInt_FMT ". The number of elements being scattered is %" PetscInt_FMT "\n", m, n));
54: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "For the Strided Scatter, fromFirst=%" PetscInt_FMT " and fromStep=%" PetscInt_FMT ".\n", fromFirst, toStep));
55: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "This produces an index (fromFirst+(n-1)*fromStep)>=m\n"));
56: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Adjust the parameterrs accordingly with -m, -n, -fromFirst, or -fromStep\n"));
57: } else {
58: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "m=%" PetscInt_FMT "\tn=%" PetscInt_FMT "\tfromFirst=%" PetscInt_FMT "\tfromStep=%" PetscInt_FMT "\ttoFirst=%" PetscInt_FMT "\ttoStep=%" PetscInt_FMT "\n", m, n, fromFirst, fromStep, toFirst, toStep));
59: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "fromFirst+(n-1)*fromStep=%" PetscInt_FMT "\ttoFirst+(n-1)*toStep=%" PetscInt_FMT "\n", fromFirst + (n - 1) * fromStep, toFirst + (n - 1) * toStep));
61: /* Build the vectors */
62: PetscCall(VecCreate(PETSC_COMM_WORLD, &Y));
63: PetscCall(VecSetSizes(Y, m, PETSC_DECIDE));
64: PetscCall(VecCreate(PETSC_COMM_WORLD, &X));
65: PetscCall(VecSetSizes(X, m, PETSC_DECIDE));
67: PetscCall(VecSetFromOptions(Y));
68: PetscCall(VecSetFromOptions(X));
69: PetscCall(VecSet(X, 2.0));
70: PetscCall(VecSet(Y, 1.0));
72: /* Build the strided index sets */
73: PetscCall(ISCreate(PETSC_COMM_WORLD, &toISStrided));
74: PetscCall(ISCreate(PETSC_COMM_WORLD, &fromISStrided));
75: PetscCall(ISSetType(toISStrided, ISSTRIDE));
76: PetscCall(ISSetType(fromISStrided, ISSTRIDE));
77: PetscCall(ISStrideSetStride(fromISStrided, n, fromFirst, fromStep));
78: PetscCall(ISStrideSetStride(toISStrided, n, toFirst, toStep));
80: /* Build the general index sets */
81: PetscCall(PetscMalloc1(n, &idx));
82: PetscCall(PetscMalloc1(n, &idy));
83: for (i = 0; i < n; i++) {
84: idx[i] = i % m;
85: idy[i] = (i + m) % m;
86: }
87: n1 = n;
88: n2 = n;
89: PetscCall(PetscSortRemoveDupsInt(&n1, idx));
90: PetscCall(PetscSortRemoveDupsInt(&n2, idy));
92: PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n1, idx, PETSC_COPY_VALUES, &toISGeneral));
93: PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n2, idy, PETSC_COPY_VALUES, &fromISGeneral));
95: /* set the mode and the insert/add parameter */
96: mode = SCATTER_FORWARD;
97: addv = ADD_VALUES;
99: /* VecScatter : Seq Strided to Seq Strided */
100: PetscCall(VecScatterCreate(X, fromISStrided, Y, toISStrided, &vscatSStoSS));
101: PetscCall(VecScatterBegin(vscatSStoSS, X, Y, addv, mode));
102: PetscCall(VecScatterEnd(vscatSStoSS, X, Y, addv, mode));
103: PetscCall(VecScatterDestroy(&vscatSStoSS));
105: /* VecScatter : Seq General to Seq Strided */
106: PetscCall(VecScatterCreate(Y, fromISGeneral, X, toISStrided, &vscatSGtoSS));
107: PetscCall(VecScatterBegin(vscatSGtoSS, Y, X, addv, mode));
108: PetscCall(VecScatterEnd(vscatSGtoSS, Y, X, addv, mode));
109: PetscCall(VecScatterDestroy(&vscatSGtoSS));
111: /* VecScatter : Seq General to Seq General */
112: PetscCall(VecScatterCreate(X, fromISGeneral, Y, toISGeneral, &vscatSGtoSG));
113: PetscCall(VecScatterBegin(vscatSGtoSG, X, Y, addv, mode));
114: PetscCall(VecScatterEnd(vscatSGtoSG, X, Y, addv, mode));
115: PetscCall(VecScatterDestroy(&vscatSGtoSG));
117: /* VecScatter : Seq Strided to Seq General */
118: PetscCall(VecScatterCreate(Y, fromISStrided, X, toISGeneral, &vscatSStoSG));
119: PetscCall(VecScatterBegin(vscatSStoSG, Y, X, addv, mode));
120: PetscCall(VecScatterEnd(vscatSStoSG, Y, X, addv, mode));
121: PetscCall(VecScatterDestroy(&vscatSStoSG));
123: /* view the results */
124: PetscCall(VecView(Y, PETSC_VIEWER_STDOUT_WORLD));
126: /* Cleanup */
127: PetscCall(VecDestroy(&X));
128: PetscCall(VecDestroy(&Y));
129: PetscCall(ISDestroy(&toISStrided));
130: PetscCall(ISDestroy(&fromISStrided));
131: PetscCall(ISDestroy(&toISGeneral));
132: PetscCall(ISDestroy(&fromISGeneral));
133: PetscCall(PetscFree(idx));
134: PetscCall(PetscFree(idy));
135: }
136: PetscCall(PetscFinalize());
137: return 0;
138: }
140: /*TEST
142: test:
143: suffix: cuda
144: args: -vec_type cuda
145: requires: cuda
147: TEST*/