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