1c4762a1bSJed Brown static char help[] = "Tests the use of MatZeroRowsColumns() for parallel matrices.\n\
2c4762a1bSJed Brown Contributed-by: Stephan Kramer <s.kramer@imperial.ac.uk>\n\n";
3c4762a1bSJed Brown
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown
main(int argc,char ** args)6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown Mat A;
9c4762a1bSJed Brown Vec x, rhs, y;
10c4762a1bSJed Brown PetscInt i, j, k, b, m = 3, n, nlocal = 2, bs = 1, Ii, J;
11c4762a1bSJed Brown PetscInt *boundary_nodes, nboundary_nodes, *boundary_indices;
12c4762a1bSJed Brown PetscMPIInt rank, size;
13c4762a1bSJed Brown PetscScalar v, v0, v1, v2, a0 = 0.1, a, rhsval, *boundary_values, diag = 1.0;
14c4762a1bSJed Brown PetscReal norm;
15c4762a1bSJed Brown char convname[64];
16c4762a1bSJed Brown PetscBool upwind = PETSC_FALSE, nonlocalBC = PETSC_FALSE, zerorhs = PETSC_TRUE, convert = PETSC_FALSE;
17c4762a1bSJed Brown
18327415f7SBarry Smith PetscFunctionBeginUser;
19c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
22c4762a1bSJed Brown n = nlocal * size;
23c4762a1bSJed Brown
249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-nonlocal_bc", &nonlocalBC, NULL));
269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL, NULL, "-diag", &diag, NULL));
279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-convname", convname, sizeof(convname), &convert));
289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-zerorhs", &zerorhs, NULL));
29c4762a1bSJed Brown
309566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
319566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m * n * bs, m * n * bs));
329566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
339566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
34c4762a1bSJed Brown
359566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, NULL, &rhs));
369566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(rhs));
379566063dSJacob Faibussowitsch PetscCall(VecSetUp(rhs));
38c4762a1bSJed Brown
39c4762a1bSJed Brown rhsval = 0.0;
40c4762a1bSJed Brown for (i = 0; i < m; i++) {
41c4762a1bSJed Brown for (j = nlocal * rank; j < nlocal * (rank + 1); j++) {
42c4762a1bSJed Brown a = a0;
43c4762a1bSJed Brown for (b = 0; b < bs; b++) {
44c4762a1bSJed Brown /* let's start with a 5-point stencil diffusion term */
459371c9d4SSatish Balay v = -1.0;
469371c9d4SSatish Balay Ii = (j + n * i) * bs + b;
479371c9d4SSatish Balay if (i > 0) {
489371c9d4SSatish Balay J = Ii - n * bs;
499371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
509371c9d4SSatish Balay }
519371c9d4SSatish Balay if (i < m - 1) {
529371c9d4SSatish Balay J = Ii + n * bs;
539371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
549371c9d4SSatish Balay }
559371c9d4SSatish Balay if (j > 0) {
569371c9d4SSatish Balay J = Ii - 1 * bs;
579371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
589371c9d4SSatish Balay }
599371c9d4SSatish Balay if (j < n - 1) {
609371c9d4SSatish Balay J = Ii + 1 * bs;
619371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
629371c9d4SSatish Balay }
639371c9d4SSatish Balay v = 4.0;
649371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, ADD_VALUES));
65c4762a1bSJed Brown if (upwind) {
66c4762a1bSJed Brown /* now add a 2nd order upwind advection term to add a little asymmetry */
67c4762a1bSJed Brown if (j > 2) {
689371c9d4SSatish Balay J = Ii - 2 * bs;
699371c9d4SSatish Balay v2 = 0.5 * a;
709371c9d4SSatish Balay v1 = -2.0 * a;
719371c9d4SSatish Balay v0 = 1.5 * a;
729566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v2, ADD_VALUES));
73c4762a1bSJed Brown } else {
74c4762a1bSJed Brown /* fall back to 1st order upwind */
759371c9d4SSatish Balay v1 = -1.0 * a;
769371c9d4SSatish Balay v0 = 1.0 * a;
77*f5c5fea7SStefano Zampini }
789371c9d4SSatish Balay if (j > 1) {
799371c9d4SSatish Balay J = Ii - 1 * bs;
809371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v1, ADD_VALUES));
819371c9d4SSatish Balay }
829566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v0, ADD_VALUES));
83c4762a1bSJed Brown a /= 10.; /* use a different velocity for the next component */
84c4762a1bSJed Brown /* add a coupling to the previous and next components */
85c4762a1bSJed Brown v = 0.5;
869371c9d4SSatish Balay if (b > 0) {
879371c9d4SSatish Balay J = Ii - 1;
889371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
899371c9d4SSatish Balay }
909371c9d4SSatish Balay if (b < bs - 1) {
919371c9d4SSatish Balay J = Ii + 1;
929371c9d4SSatish Balay PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, ADD_VALUES));
939371c9d4SSatish Balay }
94c4762a1bSJed Brown }
95c4762a1bSJed Brown /* make up some rhs */
969566063dSJacob Faibussowitsch PetscCall(VecSetValue(rhs, Ii, rhsval, INSERT_VALUES));
97c4762a1bSJed Brown rhsval += 1.0;
98c4762a1bSJed Brown }
99c4762a1bSJed Brown }
100c4762a1bSJed Brown }
1019566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1029566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
103c4762a1bSJed Brown
104c4762a1bSJed Brown if (convert) { /* Test different Mat implementations */
105c4762a1bSJed Brown Mat B;
106c4762a1bSJed Brown
1079566063dSJacob Faibussowitsch PetscCall(MatConvert(A, convname, MAT_INITIAL_MATRIX, &B));
1089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
109c4762a1bSJed Brown A = B;
110c4762a1bSJed Brown }
111c4762a1bSJed Brown
1129566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(rhs));
1139566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(rhs));
114c4762a1bSJed Brown /* set rhs to zero to simplify */
1151baa6e33SBarry Smith if (zerorhs) PetscCall(VecZeroEntries(rhs));
116c4762a1bSJed Brown
117c4762a1bSJed Brown if (nonlocalBC) {
118c4762a1bSJed Brown /*version where boundary conditions are set by processes that don't necessarily own the nodes */
119dd400576SPatrick Sanan if (rank == 0) {
120c4762a1bSJed Brown nboundary_nodes = size > m ? nlocal : m - size + nlocal;
1219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nboundary_nodes, &boundary_nodes));
122c4762a1bSJed Brown k = 0;
123ad540459SPierre Jolivet for (i = size; i < m; i++, k++) boundary_nodes[k] = n * i;
124c4762a1bSJed Brown } else if (rank < m) {
125c4762a1bSJed Brown nboundary_nodes = nlocal + 1;
1269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nboundary_nodes, &boundary_nodes));
127c4762a1bSJed Brown boundary_nodes[0] = rank * n;
128c4762a1bSJed Brown k = 1;
129c4762a1bSJed Brown } else {
130c4762a1bSJed Brown nboundary_nodes = nlocal;
1319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nboundary_nodes, &boundary_nodes));
132c4762a1bSJed Brown k = 0;
133c4762a1bSJed Brown }
134ad540459SPierre Jolivet for (j = nlocal * rank; j < nlocal * (rank + 1); j++, k++) boundary_nodes[k] = j;
135c4762a1bSJed Brown } else {
136c4762a1bSJed Brown /*version where boundary conditions are set by the node owners only */
1379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m * n, &boundary_nodes));
138c4762a1bSJed Brown k = 0;
139c4762a1bSJed Brown for (j = 0; j < n; j++) {
140c4762a1bSJed Brown Ii = j;
141c4762a1bSJed Brown if (Ii >= rank * m * nlocal && Ii < (rank + 1) * m * nlocal) boundary_nodes[k++] = Ii;
142c4762a1bSJed Brown }
143c4762a1bSJed Brown for (i = 1; i < m; i++) {
144c4762a1bSJed Brown Ii = n * i;
145c4762a1bSJed Brown if (Ii >= rank * m * nlocal && Ii < (rank + 1) * m * nlocal) boundary_nodes[k++] = Ii;
146c4762a1bSJed Brown }
147c4762a1bSJed Brown nboundary_nodes = k;
148c4762a1bSJed Brown }
149c4762a1bSJed Brown
1509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(rhs, &x));
1519566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(x));
1529566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(nboundary_nodes * bs, &boundary_indices, nboundary_nodes * bs, &boundary_values));
153c4762a1bSJed Brown for (k = 0; k < nboundary_nodes; k++) {
154c4762a1bSJed Brown Ii = boundary_nodes[k] * bs;
155c4762a1bSJed Brown v = 1.0 * boundary_nodes[k];
156c4762a1bSJed Brown for (b = 0; b < bs; b++, Ii++) {
157c4762a1bSJed Brown boundary_indices[k * bs + b] = Ii;
158c4762a1bSJed Brown boundary_values[k * bs + b] = v;
1599566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%d %" PetscInt_FMT " %f\n", rank, Ii, (double)PetscRealPart(v)));
160c4762a1bSJed Brown v += 0.1;
161c4762a1bSJed Brown }
162c4762a1bSJed Brown }
1639566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL));
1649566063dSJacob Faibussowitsch PetscCall(VecSetValues(x, nboundary_nodes * bs, boundary_indices, boundary_values, INSERT_VALUES));
1659566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x));
1669566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x));
167c4762a1bSJed Brown
168c4762a1bSJed Brown /* We can check the rhs returned by MatZeroColumns by computing y=rhs-A*x and overwriting the boundary entries with boundary values */
1699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &y));
1709566063dSJacob Faibussowitsch PetscCall(MatMult(A, x, y));
1719566063dSJacob Faibussowitsch PetscCall(VecAYPX(y, -1.0, rhs));
172c4762a1bSJed Brown for (k = 0; k < nboundary_nodes * bs; k++) boundary_values[k] *= diag;
1739566063dSJacob Faibussowitsch PetscCall(VecSetValues(y, nboundary_nodes * bs, boundary_indices, boundary_values, INSERT_VALUES));
1749566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(y));
1759566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(y));
176c4762a1bSJed Brown
1779566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "*** Matrix A and vector x:\n"));
1789566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
1799566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
180c4762a1bSJed Brown
1819566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumns(A, nboundary_nodes * bs, boundary_indices, diag, x, rhs));
1829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "*** Vector rhs returned by MatZeroRowsColumns\n"));
1839566063dSJacob Faibussowitsch PetscCall(VecView(rhs, PETSC_VIEWER_STDOUT_WORLD));
1849566063dSJacob Faibussowitsch PetscCall(VecAXPY(y, -1.0, rhs));
1859566063dSJacob Faibussowitsch PetscCall(VecNorm(y, NORM_INFINITY, &norm));
186c4762a1bSJed Brown if (norm > 1.0e-10) {
1879566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "*** Difference between rhs and y, inf-norm: %f\n", (double)norm));
1889566063dSJacob Faibussowitsch PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
189c4762a1bSJed Brown SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Bug in MatZeroRowsColumns");
190c4762a1bSJed Brown }
191c4762a1bSJed Brown
1929566063dSJacob Faibussowitsch PetscCall(PetscFree(boundary_nodes));
1939566063dSJacob Faibussowitsch PetscCall(PetscFree2(boundary_indices, boundary_values));
1949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y));
1969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rhs));
1979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
198c4762a1bSJed Brown
1999566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
200b122ec5aSJacob Faibussowitsch return 0;
201c4762a1bSJed Brown }
202c4762a1bSJed Brown
203c4762a1bSJed Brown /*TEST
204c4762a1bSJed Brown
205c4762a1bSJed Brown test:
2069fd945daSStefano Zampini diff_args: -j
207c4762a1bSJed Brown suffix: 0
208c4762a1bSJed Brown
209c4762a1bSJed Brown test:
2109fd945daSStefano Zampini diff_args: -j
211c4762a1bSJed Brown suffix: 1
212c4762a1bSJed Brown nsize: 2
213c4762a1bSJed Brown
214c4762a1bSJed Brown test:
2159fd945daSStefano Zampini diff_args: -j
216c4762a1bSJed Brown suffix: 10
217c4762a1bSJed Brown nsize: 2
218c4762a1bSJed Brown args: -bs 2 -nonlocal_bc
219c4762a1bSJed Brown
220c4762a1bSJed Brown test:
2219fd945daSStefano Zampini diff_args: -j
222c4762a1bSJed Brown suffix: 11
223c4762a1bSJed Brown nsize: 7
224c4762a1bSJed Brown args: -bs 2 -nonlocal_bc
225c4762a1bSJed Brown
226c4762a1bSJed Brown test:
2279fd945daSStefano Zampini diff_args: -j
228c4762a1bSJed Brown suffix: 12
229c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -mat_type baij
230c4762a1bSJed Brown
231c4762a1bSJed Brown test:
2329fd945daSStefano Zampini diff_args: -j
233c4762a1bSJed Brown suffix: 13
234c4762a1bSJed Brown nsize: 2
235c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -mat_type baij
236c4762a1bSJed Brown
237c4762a1bSJed Brown test:
2389fd945daSStefano Zampini diff_args: -j
239c4762a1bSJed Brown suffix: 14
240c4762a1bSJed Brown nsize: 7
241c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -mat_type baij
242c4762a1bSJed Brown
243c4762a1bSJed Brown test:
2449fd945daSStefano Zampini diff_args: -j
245c4762a1bSJed Brown suffix: 2
246c4762a1bSJed Brown nsize: 7
247c4762a1bSJed Brown
248c4762a1bSJed Brown test:
2499fd945daSStefano Zampini diff_args: -j
250c4762a1bSJed Brown suffix: 3
251c4762a1bSJed Brown args: -mat_type baij
252c4762a1bSJed Brown
253c4762a1bSJed Brown test:
2549fd945daSStefano Zampini diff_args: -j
255c4762a1bSJed Brown suffix: 4
256c4762a1bSJed Brown nsize: 2
257c4762a1bSJed Brown args: -mat_type baij
258c4762a1bSJed Brown
259c4762a1bSJed Brown test:
2609fd945daSStefano Zampini diff_args: -j
261c4762a1bSJed Brown suffix: 5
262c4762a1bSJed Brown nsize: 7
263c4762a1bSJed Brown args: -mat_type baij
264c4762a1bSJed Brown
265c4762a1bSJed Brown test:
2669fd945daSStefano Zampini diff_args: -j
267c4762a1bSJed Brown suffix: 6
268c4762a1bSJed Brown args: -bs 2 -mat_type baij
269c4762a1bSJed Brown
270c4762a1bSJed Brown test:
2719fd945daSStefano Zampini diff_args: -j
272c4762a1bSJed Brown suffix: 7
273c4762a1bSJed Brown nsize: 2
274c4762a1bSJed Brown args: -bs 2 -mat_type baij
275c4762a1bSJed Brown
276c4762a1bSJed Brown test:
2779fd945daSStefano Zampini diff_args: -j
278c4762a1bSJed Brown suffix: 8
279c4762a1bSJed Brown nsize: 7
280c4762a1bSJed Brown args: -bs 2 -mat_type baij
281c4762a1bSJed Brown
282c4762a1bSJed Brown test:
2839fd945daSStefano Zampini diff_args: -j
284c4762a1bSJed Brown suffix: 9
285c4762a1bSJed Brown args: -bs 2 -nonlocal_bc
286c4762a1bSJed Brown
287c4762a1bSJed Brown test:
2889fd945daSStefano Zampini diff_args: -j
289c4762a1bSJed Brown suffix: 15
290c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -convname shell
291c4762a1bSJed Brown
292c4762a1bSJed Brown test:
2939fd945daSStefano Zampini diff_args: -j
294c4762a1bSJed Brown suffix: 16
295c4762a1bSJed Brown nsize: 2
296c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -convname shell
297c4762a1bSJed Brown
298c4762a1bSJed Brown test:
2999fd945daSStefano Zampini diff_args: -j
300c4762a1bSJed Brown suffix: 17
301c4762a1bSJed Brown args: -bs 2 -nonlocal_bc -convname dense
302c4762a1bSJed Brown
303c4762a1bSJed Brown testset:
3049fd945daSStefano Zampini diff_args: -j
305c4762a1bSJed Brown suffix: full
306c4762a1bSJed Brown nsize: {{1 3}separate output}
307c4762a1bSJed Brown args: -diag {{0.12 -0.13}separate output} -convname {{aij shell baij}separate output} -zerorhs 0
30838a8e8c1SStefano Zampini
30938a8e8c1SStefano Zampini test:
3109fd945daSStefano Zampini diff_args: -j
31138a8e8c1SStefano Zampini requires: cuda
31238a8e8c1SStefano Zampini suffix: cusparse_1
31338a8e8c1SStefano Zampini nsize: 1
31438a8e8c1SStefano Zampini args: -diag {{0.12 -0.13}separate output} -convname {{seqaijcusparse mpiaijcusparse}separate output} -zerorhs 0 -mat_type {{seqaijcusparse mpiaijcusparse}separate output}
31538a8e8c1SStefano Zampini
31638a8e8c1SStefano Zampini test:
3179fd945daSStefano Zampini diff_args: -j
31838a8e8c1SStefano Zampini requires: cuda
31938a8e8c1SStefano Zampini suffix: cusparse_3
32038a8e8c1SStefano Zampini nsize: 3
32138a8e8c1SStefano Zampini args: -diag {{0.12 -0.13}separate output} -convname mpiaijcusparse -zerorhs 0 -mat_type mpiaijcusparse
32238a8e8c1SStefano Zampini
323c4762a1bSJed Brown TEST*/
324