xref: /petsc/src/mat/tests/ex12.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 
2 static char help[] = "Tests the use of MatZeroRows() for parallel matrices.\n\
3 This example also tests the use of MatDuplicate() for both MPIAIJ and MPIBAIJ matrices";
4 
5 #include <petscmat.h>
6 
7 extern PetscErrorCode TestMatZeroRows_Basic(Mat,IS,PetscScalar);
8 extern PetscErrorCode TestMatZeroRows_with_no_allocation(Mat,IS,PetscScalar);
9 
10 int main(int argc,char **args)
11 {
12   Mat            A;
13   PetscInt       i,j,m = 3,n,Ii,J,Imax;
14   PetscMPIInt    rank,size;
15   PetscErrorCode ierr;
16   PetscScalar    v,diag=-4.0;
17   IS             is;
18 
19   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
20   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
21   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
22   n    = 2*size;
23 
24   /* create A Square matrix for the five point stencil,YET AGAIN*/
25   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
26   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
27   CHKERRQ(MatSetFromOptions(A));
28   CHKERRQ(MatSetUp(A));
29   for (i=0; i<m; i++) {
30     for (j=2*rank; j<2*rank+2; j++) {
31       v = -1.0;  Ii = j + n*i;
32       if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
33       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
34       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
35       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
36       v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
37     }
38   }
39   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
40   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
41 
42   /* Create AN IS required by MatZeroRows() */
43   Imax = n*rank; if (Imax>= n*m -m - 1) Imax = m*n - m - 1;
44   CHKERRQ(ISCreateStride(PETSC_COMM_SELF,m,Imax,1,&is));
45 
46   CHKERRQ(TestMatZeroRows_Basic(A,is,0.0));
47   CHKERRQ(TestMatZeroRows_Basic(A,is,diag));
48 
49   CHKERRQ(TestMatZeroRows_with_no_allocation(A,is,0.0));
50   CHKERRQ(TestMatZeroRows_with_no_allocation(A,is,diag));
51 
52   CHKERRQ(MatDestroy(&A));
53 
54   /* Now Create a rectangular matrix with five point stencil (app)
55    n+size is used so that this dimension is always divisible by size.
56    This way, we can always use bs = size for any number of procs */
57   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
58   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n,m*(n+size)));
59   CHKERRQ(MatSetFromOptions(A));
60   CHKERRQ(MatSetUp(A));
61   for (i=0; i<m; i++) {
62     for (j=2*rank; j<2*rank+2; j++) {
63       v = -1.0;  Ii = j + n*i;
64       if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
65       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
66       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
67       if (j<n+size-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
68       v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
69     }
70   }
71   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
72   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
73 
74   CHKERRQ(TestMatZeroRows_Basic(A,is,0.0));
75   CHKERRQ(TestMatZeroRows_Basic(A,is,diag));
76 
77   CHKERRQ(MatDestroy(&A));
78   CHKERRQ(ISDestroy(&is));
79   ierr = PetscFinalize();
80   return ierr;
81 }
82 
83 PetscErrorCode TestMatZeroRows_Basic(Mat A,IS is,PetscScalar diag)
84 {
85   Mat            B;
86   PetscBool      keepnonzeropattern;
87 
88   /* Now copy A into B, and test it with MatZeroRows() */
89   CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&B));
90 
91   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-keep_nonzero_pattern",&keepnonzeropattern));
92   if (keepnonzeropattern) {
93     CHKERRQ(MatSetOption(B,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE));
94   }
95 
96   CHKERRQ(MatZeroRowsIS(B,is,diag,0,0));
97   CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
98   CHKERRQ(MatDestroy(&B));
99   return 0;
100 }
101 
102 PetscErrorCode TestMatZeroRows_with_no_allocation(Mat A,IS is,PetscScalar diag)
103 {
104   Mat            B;
105 
106   /* Now copy A into B, and test it with MatZeroRows() */
107   CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&B));
108   /* Set this flag after assembly. This way, it affects only MatZeroRows() */
109   CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
110 
111   CHKERRQ(MatZeroRowsIS(B,is,diag,0,0));
112   CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD));
113   CHKERRQ(MatDestroy(&B));
114   return 0;
115 }
116 
117 /*TEST
118 
119    test:
120       nsize: 2
121       filter: grep -v "MPI processes"
122 
123    test:
124       suffix: 2
125       nsize: 3
126       args: -mat_type mpibaij -mat_block_size 3
127       filter: grep -v "MPI processes"
128 
129    test:
130       suffix: 3
131       nsize: 3
132       args: -mat_type mpiaij -keep_nonzero_pattern
133       filter: grep -v "MPI processes"
134 
135    test:
136       suffix: 4
137       nsize: 3
138       args: -keep_nonzero_pattern -mat_type mpibaij -mat_block_size 3
139       filter: grep -v "MPI processes"
140 
141 TEST*/
142