xref: /petsc/src/mat/tests/ex12.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests the use of MatZeroRows() for parallel matrices.\n\
3c4762a1bSJed Brown This example also tests the use of MatDuplicate() for both MPIAIJ and MPIBAIJ matrices";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #include <petscmat.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown extern PetscErrorCode TestMatZeroRows_Basic(Mat, IS, PetscScalar);
8c4762a1bSJed Brown extern PetscErrorCode TestMatZeroRows_with_no_allocation(Mat, IS, PetscScalar);
9c4762a1bSJed Brown 
10*9371c9d4SSatish Balay int main(int argc, char **args) {
11c4762a1bSJed Brown   Mat         A;
12c4762a1bSJed Brown   PetscInt    i, j, m = 3, n, Ii, J, Imax;
13c4762a1bSJed Brown   PetscMPIInt rank, size;
14c4762a1bSJed Brown   PetscScalar v, diag = -4.0;
15c4762a1bSJed Brown   IS          is;
16c4762a1bSJed Brown 
17327415f7SBarry Smith   PetscFunctionBeginUser;
189566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
21c4762a1bSJed Brown   n = 2 * size;
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   /* create A Square matrix for the five point stencil,YET AGAIN*/
249566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
259566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
269566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
279566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
28c4762a1bSJed Brown   for (i = 0; i < m; i++) {
29c4762a1bSJed Brown     for (j = 2 * rank; j < 2 * rank + 2; j++) {
30*9371c9d4SSatish Balay       v  = -1.0;
31*9371c9d4SSatish Balay       Ii = j + n * i;
32*9371c9d4SSatish Balay       if (i > 0) {
33*9371c9d4SSatish Balay         J = Ii - n;
34*9371c9d4SSatish Balay         PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
35*9371c9d4SSatish Balay       }
36*9371c9d4SSatish Balay       if (i < m - 1) {
37*9371c9d4SSatish Balay         J = Ii + n;
38*9371c9d4SSatish Balay         PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
39*9371c9d4SSatish Balay       }
40*9371c9d4SSatish Balay       if (j > 0) {
41*9371c9d4SSatish Balay         J = Ii - 1;
42*9371c9d4SSatish Balay         PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
43*9371c9d4SSatish Balay       }
44*9371c9d4SSatish Balay       if (j < n - 1) {
45*9371c9d4SSatish Balay         J = Ii + 1;
46*9371c9d4SSatish Balay         PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
47*9371c9d4SSatish Balay       }
48*9371c9d4SSatish Balay       v = 4.0;
49*9371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
50c4762a1bSJed Brown     }
51c4762a1bSJed Brown   }
529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /* Create AN IS required by MatZeroRows() */
56*9371c9d4SSatish Balay   Imax = n * rank;
57*9371c9d4SSatish Balay   if (Imax >= n * m - m - 1) Imax = m * n - m - 1;
589566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, m, Imax, 1, &is));
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch   PetscCall(TestMatZeroRows_Basic(A, is, 0.0));
619566063dSJacob Faibussowitsch   PetscCall(TestMatZeroRows_Basic(A, is, diag));
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(TestMatZeroRows_with_no_allocation(A, is, 0.0));
649566063dSJacob Faibussowitsch   PetscCall(TestMatZeroRows_with_no_allocation(A, is, diag));
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Now Create a rectangular matrix with five point stencil (app)
69c4762a1bSJed Brown    n+size is used so that this dimension is always divisible by size.
70c4762a1bSJed Brown    This way, we can always use bs = size for any number of procs */
719566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
729566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, m * n, m * (n + size)));
739566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
749566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
75c4762a1bSJed Brown   for (i = 0; i < m; i++) {
76c4762a1bSJed Brown     for (j = 2 * rank; j < 2 * rank + 2; j++) {
77*9371c9d4SSatish Balay       v  = -1.0;
78*9371c9d4SSatish Balay       Ii = j + n * i;
79*9371c9d4SSatish Balay       if (i > 0) {
80*9371c9d4SSatish Balay         J = Ii - n;
81*9371c9d4SSatish Balay         PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
82*9371c9d4SSatish Balay       }
83*9371c9d4SSatish Balay       if (i < m - 1) {
84*9371c9d4SSatish Balay         J = Ii + n;
85*9371c9d4SSatish Balay         PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
86*9371c9d4SSatish Balay       }
87*9371c9d4SSatish Balay       if (j > 0) {
88*9371c9d4SSatish Balay         J = Ii - 1;
89*9371c9d4SSatish Balay         PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
90*9371c9d4SSatish Balay       }
91*9371c9d4SSatish Balay       if (j < n + size - 1) {
92*9371c9d4SSatish Balay         J = Ii + 1;
93*9371c9d4SSatish Balay         PetscCall(MatSetValues(A, 1, &Ii, 1, &J, &v, INSERT_VALUES));
94*9371c9d4SSatish Balay       }
95*9371c9d4SSatish Balay       v = 4.0;
96*9371c9d4SSatish Balay       PetscCall(MatSetValues(A, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
97c4762a1bSJed Brown     }
98c4762a1bSJed Brown   }
999566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1009566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
101c4762a1bSJed Brown 
1029566063dSJacob Faibussowitsch   PetscCall(TestMatZeroRows_Basic(A, is, 0.0));
1039566063dSJacob Faibussowitsch   PetscCall(TestMatZeroRows_Basic(A, is, diag));
104c4762a1bSJed Brown 
1059566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1069566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
1079566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
108b122ec5aSJacob Faibussowitsch   return 0;
109c4762a1bSJed Brown }
110c4762a1bSJed Brown 
111*9371c9d4SSatish Balay PetscErrorCode TestMatZeroRows_Basic(Mat A, IS is, PetscScalar diag) {
112c4762a1bSJed Brown   Mat       B;
113c4762a1bSJed Brown   PetscBool keepnonzeropattern;
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /* Now copy A into B, and test it with MatZeroRows() */
1169566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));
117c4762a1bSJed Brown 
1189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-keep_nonzero_pattern", &keepnonzeropattern));
1191baa6e33SBarry Smith   if (keepnonzeropattern) PetscCall(MatSetOption(B, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE));
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch   PetscCall(MatZeroRowsIS(B, is, diag, 0, 0));
1229566063dSJacob Faibussowitsch   PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD));
1239566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
124c4762a1bSJed Brown   return 0;
125c4762a1bSJed Brown }
126c4762a1bSJed Brown 
127*9371c9d4SSatish Balay PetscErrorCode TestMatZeroRows_with_no_allocation(Mat A, IS is, PetscScalar diag) {
128c4762a1bSJed Brown   Mat B;
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* Now copy A into B, and test it with MatZeroRows() */
1319566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));
132c4762a1bSJed Brown   /* Set this flag after assembly. This way, it affects only MatZeroRows() */
1339566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
134c4762a1bSJed Brown 
1359566063dSJacob Faibussowitsch   PetscCall(MatZeroRowsIS(B, is, diag, 0, 0));
1369566063dSJacob Faibussowitsch   PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD));
1379566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
138c4762a1bSJed Brown   return 0;
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown /*TEST
142c4762a1bSJed Brown 
143c4762a1bSJed Brown    test:
144c4762a1bSJed Brown       nsize: 2
1458cc725e6SPierre Jolivet       filter: grep -v " MPI process"
146c4762a1bSJed Brown 
147c4762a1bSJed Brown    test:
148c4762a1bSJed Brown       suffix: 2
149c4762a1bSJed Brown       nsize: 3
150c4762a1bSJed Brown       args: -mat_type mpibaij -mat_block_size 3
1518cc725e6SPierre Jolivet       filter: grep -v " MPI process"
152c4762a1bSJed Brown 
153c4762a1bSJed Brown    test:
154c4762a1bSJed Brown       suffix: 3
155c4762a1bSJed Brown       nsize: 3
156c4762a1bSJed Brown       args: -mat_type mpiaij -keep_nonzero_pattern
1578cc725e6SPierre Jolivet       filter: grep -v " MPI process"
158c4762a1bSJed Brown 
159c4762a1bSJed Brown    test:
160c4762a1bSJed Brown       suffix: 4
161c4762a1bSJed Brown       nsize: 3
162c4762a1bSJed Brown       args: -keep_nonzero_pattern -mat_type mpibaij -mat_block_size 3
1638cc725e6SPierre Jolivet       filter: grep -v " MPI process"
164c4762a1bSJed Brown 
165c4762a1bSJed Brown TEST*/
166