xref: /petsc/src/mat/tests/ex12.c (revision 8cc725e69398de546bdc828d7b714aa2223f5218)
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 
10c4762a1bSJed Brown int main(int argc,char **args)
11c4762a1bSJed Brown {
12c4762a1bSJed Brown   Mat            A;
13c4762a1bSJed Brown   PetscInt       i,j,m = 3,n,Ii,J,Imax;
14c4762a1bSJed Brown   PetscMPIInt    rank,size;
15c4762a1bSJed Brown   PetscScalar    v,diag=-4.0;
16c4762a1bSJed Brown   IS             is;
17c4762a1bSJed Brown 
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++) {
30c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
319566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
329566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
339566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
349566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
359566063dSJacob Faibussowitsch       v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
36c4762a1bSJed Brown     }
37c4762a1bSJed Brown   }
389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   /* Create AN IS required by MatZeroRows() */
42c4762a1bSJed Brown   Imax = n*rank; if (Imax>= n*m -m - 1) Imax = m*n - m - 1;
439566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF,m,Imax,1,&is));
44c4762a1bSJed Brown 
459566063dSJacob Faibussowitsch   PetscCall(TestMatZeroRows_Basic(A,is,0.0));
469566063dSJacob Faibussowitsch   PetscCall(TestMatZeroRows_Basic(A,is,diag));
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(TestMatZeroRows_with_no_allocation(A,is,0.0));
499566063dSJacob Faibussowitsch   PetscCall(TestMatZeroRows_with_no_allocation(A,is,diag));
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* Now Create a rectangular matrix with five point stencil (app)
54c4762a1bSJed Brown    n+size is used so that this dimension is always divisible by size.
55c4762a1bSJed Brown    This way, we can always use bs = size for any number of procs */
569566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
579566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*(n+size)));
589566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
599566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
60c4762a1bSJed Brown   for (i=0; i<m; i++) {
61c4762a1bSJed Brown     for (j=2*rank; j<2*rank+2; j++) {
62c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
639566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
649566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
659566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
669566063dSJacob Faibussowitsch       if (j<n+size-1) {J = Ii + 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
679566063dSJacob Faibussowitsch       v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
68c4762a1bSJed Brown     }
69c4762a1bSJed Brown   }
709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch   PetscCall(TestMatZeroRows_Basic(A,is,0.0));
749566063dSJacob Faibussowitsch   PetscCall(TestMatZeroRows_Basic(A,is,diag));
75c4762a1bSJed Brown 
769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
779566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
789566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
79b122ec5aSJacob Faibussowitsch   return 0;
80c4762a1bSJed Brown }
81c4762a1bSJed Brown 
82c4762a1bSJed Brown PetscErrorCode TestMatZeroRows_Basic(Mat A,IS is,PetscScalar diag)
83c4762a1bSJed Brown {
84c4762a1bSJed Brown   Mat            B;
85c4762a1bSJed Brown   PetscBool      keepnonzeropattern;
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /* Now copy A into B, and test it with MatZeroRows() */
889566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&B));
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-keep_nonzero_pattern",&keepnonzeropattern));
91c4762a1bSJed Brown   if (keepnonzeropattern) {
929566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE));
93c4762a1bSJed Brown   }
94c4762a1bSJed Brown 
959566063dSJacob Faibussowitsch   PetscCall(MatZeroRowsIS(B,is,diag,0,0));
969566063dSJacob Faibussowitsch   PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
979566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
98c4762a1bSJed Brown   return 0;
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown PetscErrorCode TestMatZeroRows_with_no_allocation(Mat A,IS is,PetscScalar diag)
102c4762a1bSJed Brown {
103c4762a1bSJed Brown   Mat            B;
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   /* Now copy A into B, and test it with MatZeroRows() */
1069566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&B));
107c4762a1bSJed Brown   /* Set this flag after assembly. This way, it affects only MatZeroRows() */
1089566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
109c4762a1bSJed Brown 
1109566063dSJacob Faibussowitsch   PetscCall(MatZeroRowsIS(B,is,diag,0,0));
1119566063dSJacob Faibussowitsch   PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
1129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
113c4762a1bSJed Brown   return 0;
114c4762a1bSJed Brown }
115c4762a1bSJed Brown 
116c4762a1bSJed Brown /*TEST
117c4762a1bSJed Brown 
118c4762a1bSJed Brown    test:
119c4762a1bSJed Brown       nsize: 2
120*8cc725e6SPierre Jolivet       filter: grep -v " MPI process"
121c4762a1bSJed Brown 
122c4762a1bSJed Brown    test:
123c4762a1bSJed Brown       suffix: 2
124c4762a1bSJed Brown       nsize: 3
125c4762a1bSJed Brown       args: -mat_type mpibaij -mat_block_size 3
126*8cc725e6SPierre Jolivet       filter: grep -v " MPI process"
127c4762a1bSJed Brown 
128c4762a1bSJed Brown    test:
129c4762a1bSJed Brown       suffix: 3
130c4762a1bSJed Brown       nsize: 3
131c4762a1bSJed Brown       args: -mat_type mpiaij -keep_nonzero_pattern
132*8cc725e6SPierre Jolivet       filter: grep -v " MPI process"
133c4762a1bSJed Brown 
134c4762a1bSJed Brown    test:
135c4762a1bSJed Brown       suffix: 4
136c4762a1bSJed Brown       nsize: 3
137c4762a1bSJed Brown       args: -keep_nonzero_pattern -mat_type mpibaij -mat_block_size 3
138*8cc725e6SPierre Jolivet       filter: grep -v " MPI process"
139c4762a1bSJed Brown 
140c4762a1bSJed Brown TEST*/
141