xref: /petsc/src/mat/tests/ex230.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static char help[] = "Example of using MatPreallocator\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown PetscErrorCode ex1_nonsquare_bs1(void)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   Mat            A,preallocator;
8c4762a1bSJed Brown   PetscInt       M,N,m,n,bs;
9c4762a1bSJed Brown 
10c4762a1bSJed Brown   /*
11c4762a1bSJed Brown      Create the Jacobian matrix
12c4762a1bSJed Brown   */
13362febeeSStefano Zampini   PetscFunctionBegin;
14c4762a1bSJed Brown   M = 10;
15c4762a1bSJed Brown   N = 8;
169566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
179566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATAIJ));
189566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
199566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(A,1));
209566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /*
23c4762a1bSJed Brown      Get the sizes of the jacobian.
24c4762a1bSJed Brown   */
259566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&m,&n));
269566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A,&bs));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /*
29c4762a1bSJed Brown      Create a preallocator matrix with sizes (local and global) matching the jacobian A.
30c4762a1bSJed Brown   */
319566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&preallocator));
329566063dSJacob Faibussowitsch   PetscCall(MatSetType(preallocator,MATPREALLOCATOR));
339566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(preallocator,m,n,M,N));
349566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(preallocator,bs));
359566063dSJacob Faibussowitsch   PetscCall(MatSetUp(preallocator));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /*
38c4762a1bSJed Brown      Insert non-zero pattern (e.g. perform a sweep over the grid).
39c4762a1bSJed Brown      You can use MatSetValues(), MatSetValuesBlocked() or MatSetValue().
40c4762a1bSJed Brown   */
41c4762a1bSJed Brown   {
42c4762a1bSJed Brown     PetscInt    ii,jj;
43c4762a1bSJed Brown     PetscScalar vv = 0.0;
44c4762a1bSJed Brown 
45c4762a1bSJed Brown     ii = 3; jj = 3;
469566063dSJacob Faibussowitsch     PetscCall(MatSetValues(preallocator,1,&ii,1,&jj,&vv,INSERT_VALUES));
47c4762a1bSJed Brown 
48c4762a1bSJed Brown     ii = 7; jj = 4;
499566063dSJacob Faibussowitsch     PetscCall(MatSetValues(preallocator,1,&ii,1,&jj,&vv,INSERT_VALUES));
50c4762a1bSJed Brown 
51c4762a1bSJed Brown     ii = 9; jj = 7;
529566063dSJacob Faibussowitsch     PetscCall(MatSetValue(preallocator,ii,jj,vv,INSERT_VALUES));
53c4762a1bSJed Brown   }
549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY));
559566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /*
58c4762a1bSJed Brown      Push the non-zero pattern defined within preallocator into A.
59c4762a1bSJed Brown      Internally this will call MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE).
60c4762a1bSJed Brown      The arg fill = PETSC_TRUE will insert zeros in the matrix A automatically.
61c4762a1bSJed Brown   */
629566063dSJacob Faibussowitsch   PetscCall(MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /*
65c4762a1bSJed Brown      We no longer require the preallocator object so destroy it.
66c4762a1bSJed Brown   */
679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&preallocator));
68c4762a1bSJed Brown 
699566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /*
72c4762a1bSJed Brown      Insert non-zero values into A.
73c4762a1bSJed Brown   */
74c4762a1bSJed Brown   {
75c4762a1bSJed Brown     PetscInt    ii,jj;
76c4762a1bSJed Brown     PetscScalar vv;
77c4762a1bSJed Brown 
78c4762a1bSJed Brown     ii = 3; jj = 3; vv = 0.3;
799566063dSJacob Faibussowitsch     PetscCall(MatSetValue(A,ii,jj,vv,INSERT_VALUES));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown     ii = 7; jj = 4; vv = 3.3;
829566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,1,&ii,1,&jj,&vv,INSERT_VALUES));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown     ii = 9; jj = 7; vv = 4.3;
859566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A,1,&ii,1,&jj,&vv,INSERT_VALUES));
86c4762a1bSJed Brown   }
879566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
889566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
91c4762a1bSJed Brown 
929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
93c4762a1bSJed Brown   PetscFunctionReturn(0);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown PetscErrorCode ex2_square_bsvariable(void)
97c4762a1bSJed Brown {
98c4762a1bSJed Brown   Mat            A,preallocator;
99c4762a1bSJed Brown   PetscInt       M,N,m,n,bs = 1;
100c4762a1bSJed Brown 
101362febeeSStefano Zampini   PetscFunctionBegin;
1029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-block_size",&bs,NULL));
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   /*
105c4762a1bSJed Brown      Create the Jacobian matrix.
106c4762a1bSJed Brown   */
107c4762a1bSJed Brown   M = 10 * bs;
108c4762a1bSJed Brown   N = 10 * bs;
1099566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
1109566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATAIJ));
1119566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
1129566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(A,bs));
1139566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /*
116c4762a1bSJed Brown      Get the sizes of the jacobian.
117c4762a1bSJed Brown   */
1189566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&m,&n));
1199566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A,&bs));
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /*
122c4762a1bSJed Brown      Create a preallocator matrix with dimensions matching the jacobian A.
123c4762a1bSJed Brown   */
1249566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&preallocator));
1259566063dSJacob Faibussowitsch   PetscCall(MatSetType(preallocator,MATPREALLOCATOR));
1269566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(preallocator,m,n,M,N));
1279566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(preallocator,bs));
1289566063dSJacob Faibussowitsch   PetscCall(MatSetUp(preallocator));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /*
131c4762a1bSJed Brown      Insert non-zero pattern (e.g. perform a sweep over the grid).
132c4762a1bSJed Brown      You can use MatSetValues(), MatSetValuesBlocked() or MatSetValue().
133c4762a1bSJed Brown   */
134c4762a1bSJed Brown   {
135c4762a1bSJed Brown     PetscInt  ii,jj;
136c4762a1bSJed Brown     PetscScalar *vv;
137c4762a1bSJed Brown 
1389566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(bs*bs,&vv));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown     ii = 0; jj = 9;
1419566063dSJacob Faibussowitsch     PetscCall(MatSetValue(preallocator,ii,jj,vv[0],INSERT_VALUES));
142c4762a1bSJed Brown 
143c4762a1bSJed Brown     ii = 0; jj = 0;
1449566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown     ii = 2; jj = 4;
1479566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES));
148c4762a1bSJed Brown 
149c4762a1bSJed Brown     ii = 4; jj = 2;
1509566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown     ii = 4; jj = 4;
1539566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES));
154c4762a1bSJed Brown 
155c4762a1bSJed Brown     ii = 9; jj = 9;
1569566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES));
157c4762a1bSJed Brown 
1589566063dSJacob Faibussowitsch     PetscCall(PetscFree(vv));
159c4762a1bSJed Brown   }
1609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY));
1619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /*
164c4762a1bSJed Brown      Push non-zero pattern defined from preallocator into A.
165c4762a1bSJed Brown      Internally this will call MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE).
166c4762a1bSJed Brown      The arg fill = PETSC_TRUE will insert zeros in the matrix A automatically.
167c4762a1bSJed Brown   */
1689566063dSJacob Faibussowitsch   PetscCall(MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A));
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /*
171c4762a1bSJed Brown      We no longer require the preallocator object so destroy it.
172c4762a1bSJed Brown   */
1739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&preallocator));
174c4762a1bSJed Brown 
1759566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   {
178c4762a1bSJed Brown     PetscInt    ii,jj;
179c4762a1bSJed Brown     PetscScalar *vv;
180c4762a1bSJed Brown 
1819566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(bs*bs,&vv));
182c4762a1bSJed Brown     for (ii=0; ii<bs*bs; ii++) {
183c4762a1bSJed Brown       vv[ii] = (PetscReal)(ii+1);
184c4762a1bSJed Brown     }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown     ii = 0; jj = 9;
1879566063dSJacob Faibussowitsch     PetscCall(MatSetValue(A,ii,jj,33.3,INSERT_VALUES));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown     ii = 0; jj = 0;
1909566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES));
191c4762a1bSJed Brown 
192c4762a1bSJed Brown     ii = 2; jj = 4;
1939566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES));
194c4762a1bSJed Brown 
195c4762a1bSJed Brown     ii = 4; jj = 2;
1969566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES));
197c4762a1bSJed Brown 
198c4762a1bSJed Brown     ii = 4; jj = 4;
1999566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES));
200c4762a1bSJed Brown 
201c4762a1bSJed Brown     ii = 9; jj = 9;
2029566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES));
203c4762a1bSJed Brown 
2049566063dSJacob Faibussowitsch     PetscCall(PetscFree(vv));
205c4762a1bSJed Brown   }
2069566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
2079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
208c4762a1bSJed Brown 
2099566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
210c4762a1bSJed Brown 
2119566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
212c4762a1bSJed Brown   PetscFunctionReturn(0);
213c4762a1bSJed Brown }
214c4762a1bSJed Brown 
215c4762a1bSJed Brown int main(int argc, char **args)
216c4762a1bSJed Brown {
217c4762a1bSJed Brown   PetscInt testid = 0;
218*327415f7SBarry Smith   PetscFunctionBeginUser;
2199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
2209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-test_id",&testid,NULL));
221c4762a1bSJed Brown   switch (testid) {
222c4762a1bSJed Brown     case 0:
2239566063dSJacob Faibussowitsch       PetscCall(ex1_nonsquare_bs1());
224c4762a1bSJed Brown       break;
225c4762a1bSJed Brown     case 1:
2269566063dSJacob Faibussowitsch       PetscCall(ex2_square_bsvariable());
227c4762a1bSJed Brown       break;
228c4762a1bSJed Brown     default:
229c4762a1bSJed Brown       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Invalid value for -test_id. Must be {0,1}");
230c4762a1bSJed Brown   }
2319566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
232b122ec5aSJacob Faibussowitsch   return 0;
233c4762a1bSJed Brown }
234c4762a1bSJed Brown 
235c4762a1bSJed Brown /*TEST
236c4762a1bSJed Brown 
237c4762a1bSJed Brown    test:
238c4762a1bSJed Brown      suffix: t0_a_aij
239c4762a1bSJed Brown      nsize: 1
240c4762a1bSJed Brown      args: -test_id 0 -mat_type aij
241c4762a1bSJed Brown 
242c4762a1bSJed Brown    test:
243c4762a1bSJed Brown      suffix: t0_b_aij
244c4762a1bSJed Brown      nsize: 6
245c4762a1bSJed Brown      args: -test_id 0 -mat_type aij
246c4762a1bSJed Brown 
247c4762a1bSJed Brown    test:
248c4762a1bSJed Brown      suffix: t1_a_aij
249c4762a1bSJed Brown      nsize: 1
250c4762a1bSJed Brown      args: -test_id 1 -mat_type aij
251c4762a1bSJed Brown 
252c4762a1bSJed Brown    test:
253c4762a1bSJed Brown      suffix: t2_a_baij
254c4762a1bSJed Brown      nsize: 1
255c4762a1bSJed Brown      args: -test_id 1 -mat_type baij
256c4762a1bSJed Brown 
257c4762a1bSJed Brown    test:
258c4762a1bSJed Brown      suffix: t3_a_sbaij
259c4762a1bSJed Brown      nsize: 1
260c4762a1bSJed Brown      args: -test_id 1 -mat_type sbaij
261c4762a1bSJed Brown 
262c4762a1bSJed Brown    test:
263c4762a1bSJed Brown      suffix: t4_a_aij_bs3
264c4762a1bSJed Brown      nsize: 1
265c4762a1bSJed Brown      args: -test_id 1 -mat_type aij -block_size 3
266c4762a1bSJed Brown 
267c4762a1bSJed Brown    test:
268c4762a1bSJed Brown      suffix: t5_a_baij_bs3
269c4762a1bSJed Brown      nsize: 1
270c4762a1bSJed Brown      args: -test_id 1 -mat_type baij -block_size 3
271c4762a1bSJed Brown 
272c4762a1bSJed Brown    test:
273c4762a1bSJed Brown      suffix: t6_a_sbaij_bs3
274c4762a1bSJed Brown      nsize: 1
275c4762a1bSJed Brown      args: -test_id 1 -mat_type sbaij -block_size 3
276c4762a1bSJed Brown 
277c4762a1bSJed Brown    test:
278c4762a1bSJed Brown      suffix: t4_b_aij_bs3
279c4762a1bSJed Brown      nsize: 6
280c4762a1bSJed Brown      args: -test_id 1 -mat_type aij -block_size 3
281c4762a1bSJed Brown 
282c4762a1bSJed Brown    test:
283c4762a1bSJed Brown      suffix: t5_b_baij_bs3
284c4762a1bSJed Brown      nsize: 6
285c4762a1bSJed Brown      args: -test_id 1 -mat_type baij -block_size 3
286c4762a1bSJed Brown      filter: grep -v Mat_
287c4762a1bSJed Brown 
288c4762a1bSJed Brown    test:
289c4762a1bSJed Brown      suffix: t6_b_sbaij_bs3
290c4762a1bSJed Brown      nsize: 6
291c4762a1bSJed Brown      args: -test_id 1 -mat_type sbaij -block_size 3
292c4762a1bSJed Brown      filter: grep -v Mat_
293c4762a1bSJed Brown 
294c4762a1bSJed Brown TEST*/
295