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