xref: /petsc/src/mat/tests/ex18.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Tests the use of MatZeroRowsColumns() for parallel matrices.\n\
2*c4762a1bSJed Brown Contributed-by: Stephan Kramer <s.kramer@imperial.ac.uk>\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown int main(int argc,char **args)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   Mat            A;
9*c4762a1bSJed Brown   Vec            x, rhs, y;
10*c4762a1bSJed Brown   PetscInt       i,j,k,b,m = 3,n,nlocal=2,bs=1,Ii,J;
11*c4762a1bSJed Brown   PetscInt       *boundary_nodes, nboundary_nodes, *boundary_indices;
12*c4762a1bSJed Brown   PetscMPIInt    rank,size;
13*c4762a1bSJed Brown   PetscErrorCode ierr;
14*c4762a1bSJed Brown   PetscScalar    v,v0,v1,v2,a0=0.1,a,rhsval, *boundary_values,diag = 1.0;
15*c4762a1bSJed Brown   PetscReal      norm;
16*c4762a1bSJed Brown   char           convname[64];
17*c4762a1bSJed Brown   PetscBool      upwind = PETSC_FALSE, nonlocalBC = PETSC_FALSE, zerorhs = PETSC_TRUE, convert = PETSC_FALSE;
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 = nlocal*size;
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL, "-bs", &bs, NULL);CHKERRQ(ierr);
25*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL, "-nonlocal_bc", &nonlocalBC, NULL);CHKERRQ(ierr);
26*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL, "-diag", &diag, NULL);CHKERRQ(ierr);
27*c4762a1bSJed Brown   ierr = PetscOptionsGetString(NULL,NULL,"-convname",convname,sizeof(convname),&convert);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL, "-zerorhs", &zerorhs, NULL);CHKERRQ(ierr);
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m*n*bs,m*n*bs);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD, &rhs);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = VecSetSizes(rhs, PETSC_DECIDE, m*n*bs);CHKERRQ(ierr);
37*c4762a1bSJed Brown   ierr = VecSetFromOptions(rhs);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = VecSetUp(rhs);CHKERRQ(ierr);
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown   rhsval = 0.0;
41*c4762a1bSJed Brown   for (i=0; i<m; i++) {
42*c4762a1bSJed Brown     for (j=nlocal*rank; j<nlocal*(rank+1); j++) {
43*c4762a1bSJed Brown       a = a0;
44*c4762a1bSJed Brown       for (b=0; b<bs; b++) {
45*c4762a1bSJed Brown         /* let's start with a 5-point stencil diffusion term */
46*c4762a1bSJed Brown         v = -1.0;  Ii = (j + n*i)*bs + b;
47*c4762a1bSJed Brown         if (i>0)   {J = Ii - n*bs; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);}
48*c4762a1bSJed Brown         if (i<m-1) {J = Ii + n*bs; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);}
49*c4762a1bSJed Brown         if (j>0)   {J = Ii - 1*bs; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);}
50*c4762a1bSJed Brown         if (j<n-1) {J = Ii + 1*bs; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);}
51*c4762a1bSJed Brown         v = 4.0; ierr = MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES);CHKERRQ(ierr);
52*c4762a1bSJed Brown         if (upwind) {
53*c4762a1bSJed Brown           /* now add a 2nd order upwind advection term to add a little asymmetry */
54*c4762a1bSJed Brown           if (j>2) {
55*c4762a1bSJed Brown             J = Ii-2*bs; v2 = 0.5*a; v1 = -2.0*a; v0 = 1.5*a;
56*c4762a1bSJed Brown             ierr = MatSetValues(A,1,&Ii,1,&J,&v2,ADD_VALUES);CHKERRQ(ierr);
57*c4762a1bSJed Brown           } else {
58*c4762a1bSJed Brown             /* fall back to 1st order upwind */
59*c4762a1bSJed Brown             v1 = -1.0*a; v0 = 1.0*a;
60*c4762a1bSJed Brown           };
61*c4762a1bSJed Brown           if (j>1) {J = Ii-1*bs; ierr = MatSetValues(A,1,&Ii,1,&J,&v1,ADD_VALUES);CHKERRQ(ierr);}
62*c4762a1bSJed Brown           ierr = MatSetValues(A,1,&Ii,1,&Ii,&v0,ADD_VALUES);CHKERRQ(ierr);
63*c4762a1bSJed Brown           a /= 10.; /* use a different velocity for the next component */
64*c4762a1bSJed Brown           /* add a coupling to the previous and next components */
65*c4762a1bSJed Brown           v = 0.5;
66*c4762a1bSJed Brown           if (b>0) {J = Ii - 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);}
67*c4762a1bSJed Brown           if (b<bs-1) {J = Ii + 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);}
68*c4762a1bSJed Brown         }
69*c4762a1bSJed Brown         /* make up some rhs */
70*c4762a1bSJed Brown         ierr = VecSetValue(rhs, Ii, rhsval, INSERT_VALUES);CHKERRQ(ierr);
71*c4762a1bSJed Brown         rhsval += 1.0;
72*c4762a1bSJed Brown       }
73*c4762a1bSJed Brown     }
74*c4762a1bSJed Brown   }
75*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown   if (convert) { /* Test different Mat implementations */
79*c4762a1bSJed Brown     Mat B;
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown     ierr = MatConvert(A,convname,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
82*c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
83*c4762a1bSJed Brown     A    = B;
84*c4762a1bSJed Brown   }
85*c4762a1bSJed Brown 
86*c4762a1bSJed Brown   ierr = VecAssemblyBegin(rhs);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = VecAssemblyEnd(rhs);CHKERRQ(ierr);
88*c4762a1bSJed Brown   /* set rhs to zero to simplify */
89*c4762a1bSJed Brown   if (zerorhs) {
90*c4762a1bSJed Brown     ierr = VecZeroEntries(rhs);CHKERRQ(ierr);
91*c4762a1bSJed Brown   }
92*c4762a1bSJed Brown 
93*c4762a1bSJed Brown   if (nonlocalBC) {
94*c4762a1bSJed Brown     /*version where boundary conditions are set by processes that don't necessarily own the nodes */
95*c4762a1bSJed Brown     if (!rank) {
96*c4762a1bSJed Brown       nboundary_nodes = size>m ? nlocal : m-size+nlocal;
97*c4762a1bSJed Brown       ierr = PetscMalloc1(nboundary_nodes,&boundary_nodes);CHKERRQ(ierr);
98*c4762a1bSJed Brown       k = 0;
99*c4762a1bSJed Brown       for (i=size; i<m; i++,k++) {boundary_nodes[k] = n*i;};
100*c4762a1bSJed Brown     } else if (rank < m) {
101*c4762a1bSJed Brown       nboundary_nodes = nlocal+1;
102*c4762a1bSJed Brown       ierr = PetscMalloc1(nboundary_nodes,&boundary_nodes);CHKERRQ(ierr);
103*c4762a1bSJed Brown       boundary_nodes[0] = rank*n;
104*c4762a1bSJed Brown       k = 1;
105*c4762a1bSJed Brown     } else {
106*c4762a1bSJed Brown       nboundary_nodes = nlocal;
107*c4762a1bSJed Brown       ierr = PetscMalloc1(nboundary_nodes,&boundary_nodes);CHKERRQ(ierr);
108*c4762a1bSJed Brown       k = 0;
109*c4762a1bSJed Brown     }
110*c4762a1bSJed Brown     for (j=nlocal*rank; j<nlocal*(rank+1); j++,k++) {boundary_nodes[k] = j;};
111*c4762a1bSJed Brown   } else {
112*c4762a1bSJed Brown     /*version where boundary conditions are set by the node owners only */
113*c4762a1bSJed Brown     ierr = PetscMalloc1(m*n,&boundary_nodes);CHKERRQ(ierr);
114*c4762a1bSJed Brown     k=0;
115*c4762a1bSJed Brown     for (j=0; j<n; j++) {
116*c4762a1bSJed Brown       Ii = j;
117*c4762a1bSJed Brown       if (Ii>=rank*m*nlocal && Ii<(rank+1)*m*nlocal) boundary_nodes[k++] = Ii;
118*c4762a1bSJed Brown     }
119*c4762a1bSJed Brown     for (i=1; i<m; i++) {
120*c4762a1bSJed Brown       Ii = n*i;
121*c4762a1bSJed Brown       if (Ii>=rank*m*nlocal && Ii<(rank+1)*m*nlocal) boundary_nodes[k++] = Ii;
122*c4762a1bSJed Brown     }
123*c4762a1bSJed Brown     nboundary_nodes = k;
124*c4762a1bSJed Brown   }
125*c4762a1bSJed Brown 
126*c4762a1bSJed Brown   ierr = VecDuplicate(rhs, &x);CHKERRQ(ierr);
127*c4762a1bSJed Brown   ierr = VecZeroEntries(x);CHKERRQ(ierr);
128*c4762a1bSJed Brown   ierr = PetscMalloc2(nboundary_nodes*bs,&boundary_indices,nboundary_nodes*bs,&boundary_values);CHKERRQ(ierr);
129*c4762a1bSJed Brown   for (k=0; k<nboundary_nodes; k++) {
130*c4762a1bSJed Brown     Ii = boundary_nodes[k]*bs;
131*c4762a1bSJed Brown     v = 1.0*boundary_nodes[k];
132*c4762a1bSJed Brown     for (b=0; b<bs; b++, Ii++) {
133*c4762a1bSJed Brown       boundary_indices[k*bs+b] = Ii;
134*c4762a1bSJed Brown       boundary_values[k*bs+b] = v;
135*c4762a1bSJed Brown       ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD, "%d %D %f\n", rank, Ii, (double)PetscRealPart(v));CHKERRQ(ierr);
136*c4762a1bSJed Brown       v += 0.1;
137*c4762a1bSJed Brown     }
138*c4762a1bSJed Brown   }
139*c4762a1bSJed Brown   ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD, NULL);CHKERRQ(ierr);
140*c4762a1bSJed Brown   ierr = VecSetValues(x, nboundary_nodes*bs, boundary_indices, boundary_values, INSERT_VALUES);CHKERRQ(ierr);
141*c4762a1bSJed Brown   ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
143*c4762a1bSJed Brown 
144*c4762a1bSJed Brown   /* We can check the rhs returned by MatZeroColumns by computing y=rhs-A*x  and overwriting the boundary entries with boundary values */
145*c4762a1bSJed Brown   ierr = VecDuplicate(x, &y);CHKERRQ(ierr);
146*c4762a1bSJed Brown   ierr = MatMult(A, x, y);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = VecAYPX(y, -1.0, rhs);CHKERRQ(ierr);
148*c4762a1bSJed Brown   for (k=0; k<nboundary_nodes*bs; k++) boundary_values[k] *= diag;
149*c4762a1bSJed Brown   ierr = VecSetValues(y, nboundary_nodes*bs, boundary_indices, boundary_values, INSERT_VALUES);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = VecAssemblyBegin(y);CHKERRQ(ierr);
151*c4762a1bSJed Brown   ierr = VecAssemblyEnd(y);CHKERRQ(ierr);
152*c4762a1bSJed Brown 
153*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "*** Matrix A and vector x:\n");CHKERRQ(ierr);
154*c4762a1bSJed Brown   ierr = MatView(A, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
155*c4762a1bSJed Brown   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
156*c4762a1bSJed Brown 
157*c4762a1bSJed Brown   ierr = MatZeroRowsColumns(A, nboundary_nodes*bs, boundary_indices, diag, x, rhs);CHKERRQ(ierr);
158*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "*** Vector rhs returned by MatZeroRowsColumns\n");CHKERRQ(ierr);
159*c4762a1bSJed Brown   ierr = VecView(rhs,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
160*c4762a1bSJed Brown   ierr = VecAXPY(y, -1.0, rhs);CHKERRQ(ierr);
161*c4762a1bSJed Brown   ierr = VecNorm(y, NORM_INFINITY, &norm);CHKERRQ(ierr);
162*c4762a1bSJed Brown   if (norm > 1.0e-10) {
163*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "*** Difference between rhs and y, inf-norm: %f\n", (double)norm);CHKERRQ(ierr);
164*c4762a1bSJed Brown     ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
165*c4762a1bSJed Brown     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Bug in MatZeroRowsColumns");
166*c4762a1bSJed Brown   }
167*c4762a1bSJed Brown 
168*c4762a1bSJed Brown   ierr = PetscFree(boundary_nodes);CHKERRQ(ierr);
169*c4762a1bSJed Brown   ierr = PetscFree2(boundary_indices,boundary_values);CHKERRQ(ierr);
170*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
171*c4762a1bSJed Brown   ierr = VecDestroy(&y);CHKERRQ(ierr);
172*c4762a1bSJed Brown   ierr = VecDestroy(&rhs);CHKERRQ(ierr);
173*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
174*c4762a1bSJed Brown 
175*c4762a1bSJed Brown   ierr = PetscFinalize();
176*c4762a1bSJed Brown   return ierr;
177*c4762a1bSJed Brown }
178*c4762a1bSJed Brown 
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown /*TEST
181*c4762a1bSJed Brown 
182*c4762a1bSJed Brown    test:
183*c4762a1bSJed Brown       suffix: 0
184*c4762a1bSJed Brown 
185*c4762a1bSJed Brown    test:
186*c4762a1bSJed Brown       suffix: 1
187*c4762a1bSJed Brown       nsize: 2
188*c4762a1bSJed Brown 
189*c4762a1bSJed Brown    test:
190*c4762a1bSJed Brown       suffix: 10
191*c4762a1bSJed Brown       nsize: 2
192*c4762a1bSJed Brown       args: -bs 2 -nonlocal_bc
193*c4762a1bSJed Brown 
194*c4762a1bSJed Brown    test:
195*c4762a1bSJed Brown       suffix: 11
196*c4762a1bSJed Brown       nsize: 7
197*c4762a1bSJed Brown       args: -bs 2 -nonlocal_bc
198*c4762a1bSJed Brown 
199*c4762a1bSJed Brown    test:
200*c4762a1bSJed Brown       suffix: 12
201*c4762a1bSJed Brown       args: -bs 2 -nonlocal_bc -mat_type baij
202*c4762a1bSJed Brown 
203*c4762a1bSJed Brown    test:
204*c4762a1bSJed Brown       suffix: 13
205*c4762a1bSJed Brown       nsize: 2
206*c4762a1bSJed Brown       args: -bs 2 -nonlocal_bc -mat_type baij
207*c4762a1bSJed Brown 
208*c4762a1bSJed Brown    test:
209*c4762a1bSJed Brown       suffix: 14
210*c4762a1bSJed Brown       nsize: 7
211*c4762a1bSJed Brown       args: -bs 2 -nonlocal_bc -mat_type baij
212*c4762a1bSJed Brown 
213*c4762a1bSJed Brown    test:
214*c4762a1bSJed Brown       suffix: 2
215*c4762a1bSJed Brown       nsize: 7
216*c4762a1bSJed Brown 
217*c4762a1bSJed Brown    test:
218*c4762a1bSJed Brown       suffix: 3
219*c4762a1bSJed Brown       args: -mat_type baij
220*c4762a1bSJed Brown 
221*c4762a1bSJed Brown    test:
222*c4762a1bSJed Brown       suffix: 4
223*c4762a1bSJed Brown       nsize: 2
224*c4762a1bSJed Brown       args: -mat_type baij
225*c4762a1bSJed Brown 
226*c4762a1bSJed Brown    test:
227*c4762a1bSJed Brown       suffix: 5
228*c4762a1bSJed Brown       nsize: 7
229*c4762a1bSJed Brown       args: -mat_type baij
230*c4762a1bSJed Brown 
231*c4762a1bSJed Brown    test:
232*c4762a1bSJed Brown       suffix: 6
233*c4762a1bSJed Brown       args: -bs 2 -mat_type baij
234*c4762a1bSJed Brown 
235*c4762a1bSJed Brown    test:
236*c4762a1bSJed Brown       suffix: 7
237*c4762a1bSJed Brown       nsize: 2
238*c4762a1bSJed Brown       args: -bs 2 -mat_type baij
239*c4762a1bSJed Brown 
240*c4762a1bSJed Brown    test:
241*c4762a1bSJed Brown       suffix: 8
242*c4762a1bSJed Brown       nsize: 7
243*c4762a1bSJed Brown       args: -bs 2 -mat_type baij
244*c4762a1bSJed Brown 
245*c4762a1bSJed Brown    test:
246*c4762a1bSJed Brown       suffix: 9
247*c4762a1bSJed Brown       args: -bs 2 -nonlocal_bc
248*c4762a1bSJed Brown 
249*c4762a1bSJed Brown    test:
250*c4762a1bSJed Brown       suffix: 15
251*c4762a1bSJed Brown       args: -bs 2 -nonlocal_bc -convname shell
252*c4762a1bSJed Brown 
253*c4762a1bSJed Brown    test:
254*c4762a1bSJed Brown       suffix: 16
255*c4762a1bSJed Brown       nsize: 2
256*c4762a1bSJed Brown       args: -bs 2 -nonlocal_bc -convname shell
257*c4762a1bSJed Brown 
258*c4762a1bSJed Brown    test:
259*c4762a1bSJed Brown       suffix: 17
260*c4762a1bSJed Brown       args: -bs 2 -nonlocal_bc -convname dense
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown    testset:
263*c4762a1bSJed Brown       suffix: full
264*c4762a1bSJed Brown       nsize: {{1 3}separate output}
265*c4762a1bSJed Brown       args: -diag {{0.12 -0.13}separate output} -convname {{aij shell baij}separate output} -zerorhs 0
266*c4762a1bSJed Brown TEST*/
267