xref: /petsc/src/mat/tests/ex230.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Example of using MatPreallocator\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5*9371c9d4SSatish Balay PetscErrorCode ex1_nonsquare_bs1(void) {
6c4762a1bSJed Brown   Mat      A, preallocator;
7c4762a1bSJed Brown   PetscInt M, N, m, n, bs;
8c4762a1bSJed Brown 
9c4762a1bSJed Brown   /*
10c4762a1bSJed Brown      Create the Jacobian matrix
11c4762a1bSJed Brown   */
12362febeeSStefano Zampini   PetscFunctionBegin;
13c4762a1bSJed Brown   M = 10;
14c4762a1bSJed Brown   N = 8;
159566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
169566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
179566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N));
189566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(A, 1));
199566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   /*
22c4762a1bSJed Brown      Get the sizes of the jacobian.
23c4762a1bSJed Brown   */
249566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
259566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   /*
28c4762a1bSJed Brown      Create a preallocator matrix with sizes (local and global) matching the jacobian A.
29c4762a1bSJed Brown   */
309566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator));
319566063dSJacob Faibussowitsch   PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
329566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(preallocator, m, n, M, N));
339566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(preallocator, bs));
349566063dSJacob Faibussowitsch   PetscCall(MatSetUp(preallocator));
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   /*
37c4762a1bSJed Brown      Insert non-zero pattern (e.g. perform a sweep over the grid).
38c4762a1bSJed Brown      You can use MatSetValues(), MatSetValuesBlocked() or MatSetValue().
39c4762a1bSJed Brown   */
40c4762a1bSJed Brown   {
41c4762a1bSJed Brown     PetscInt    ii, jj;
42c4762a1bSJed Brown     PetscScalar vv = 0.0;
43c4762a1bSJed Brown 
44*9371c9d4SSatish Balay     ii = 3;
45*9371c9d4SSatish Balay     jj = 3;
469566063dSJacob Faibussowitsch     PetscCall(MatSetValues(preallocator, 1, &ii, 1, &jj, &vv, INSERT_VALUES));
47c4762a1bSJed Brown 
48*9371c9d4SSatish Balay     ii = 7;
49*9371c9d4SSatish Balay     jj = 4;
509566063dSJacob Faibussowitsch     PetscCall(MatSetValues(preallocator, 1, &ii, 1, &jj, &vv, INSERT_VALUES));
51c4762a1bSJed Brown 
52*9371c9d4SSatish Balay     ii = 9;
53*9371c9d4SSatish Balay     jj = 7;
549566063dSJacob Faibussowitsch     PetscCall(MatSetValue(preallocator, ii, jj, vv, INSERT_VALUES));
55c4762a1bSJed Brown   }
569566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /*
60c4762a1bSJed Brown      Push the non-zero pattern defined within preallocator into A.
61c4762a1bSJed Brown      Internally this will call MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE).
62c4762a1bSJed Brown      The arg fill = PETSC_TRUE will insert zeros in the matrix A automatically.
63c4762a1bSJed Brown   */
649566063dSJacob Faibussowitsch   PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, A));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /*
67c4762a1bSJed Brown      We no longer require the preallocator object so destroy it.
68c4762a1bSJed Brown   */
699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&preallocator));
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /*
74c4762a1bSJed Brown      Insert non-zero values into A.
75c4762a1bSJed Brown   */
76c4762a1bSJed Brown   {
77c4762a1bSJed Brown     PetscInt    ii, jj;
78c4762a1bSJed Brown     PetscScalar vv;
79c4762a1bSJed Brown 
80*9371c9d4SSatish Balay     ii = 3;
81*9371c9d4SSatish Balay     jj = 3;
82*9371c9d4SSatish Balay     vv = 0.3;
839566063dSJacob Faibussowitsch     PetscCall(MatSetValue(A, ii, jj, vv, INSERT_VALUES));
84c4762a1bSJed Brown 
85*9371c9d4SSatish Balay     ii = 7;
86*9371c9d4SSatish Balay     jj = 4;
87*9371c9d4SSatish Balay     vv = 3.3;
889566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &ii, 1, &jj, &vv, INSERT_VALUES));
89c4762a1bSJed Brown 
90*9371c9d4SSatish Balay     ii = 9;
91*9371c9d4SSatish Balay     jj = 7;
92*9371c9d4SSatish Balay     vv = 4.3;
939566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &ii, 1, &jj, &vv, INSERT_VALUES));
94c4762a1bSJed Brown   }
959566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
97c4762a1bSJed Brown 
989566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
99c4762a1bSJed Brown 
1009566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
101c4762a1bSJed Brown   PetscFunctionReturn(0);
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104*9371c9d4SSatish Balay PetscErrorCode ex2_square_bsvariable(void) {
105c4762a1bSJed Brown   Mat      A, preallocator;
106c4762a1bSJed Brown   PetscInt M, N, m, n, bs = 1;
107c4762a1bSJed Brown 
108362febeeSStefano Zampini   PetscFunctionBegin;
1099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-block_size", &bs, NULL));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /*
112c4762a1bSJed Brown      Create the Jacobian matrix.
113c4762a1bSJed Brown   */
114c4762a1bSJed Brown   M = 10 * bs;
115c4762a1bSJed Brown   N = 10 * bs;
1169566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1179566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
1189566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N));
1199566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(A, bs));
1209566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   /*
123c4762a1bSJed Brown      Get the sizes of the jacobian.
124c4762a1bSJed Brown   */
1259566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &m, &n));
1269566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   /*
129c4762a1bSJed Brown      Create a preallocator matrix with dimensions matching the jacobian A.
130c4762a1bSJed Brown   */
1319566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator));
1329566063dSJacob Faibussowitsch   PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
1339566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(preallocator, m, n, M, N));
1349566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(preallocator, bs));
1359566063dSJacob Faibussowitsch   PetscCall(MatSetUp(preallocator));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /*
138c4762a1bSJed Brown      Insert non-zero pattern (e.g. perform a sweep over the grid).
139c4762a1bSJed Brown      You can use MatSetValues(), MatSetValuesBlocked() or MatSetValue().
140c4762a1bSJed Brown   */
141c4762a1bSJed Brown   {
142c4762a1bSJed Brown     PetscInt     ii, jj;
143c4762a1bSJed Brown     PetscScalar *vv;
144c4762a1bSJed Brown 
1459566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(bs * bs, &vv));
146c4762a1bSJed Brown 
147*9371c9d4SSatish Balay     ii = 0;
148*9371c9d4SSatish Balay     jj = 9;
1499566063dSJacob Faibussowitsch     PetscCall(MatSetValue(preallocator, ii, jj, vv[0], INSERT_VALUES));
150c4762a1bSJed Brown 
151*9371c9d4SSatish Balay     ii = 0;
152*9371c9d4SSatish Balay     jj = 0;
1539566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(preallocator, 1, &ii, 1, &jj, vv, INSERT_VALUES));
154c4762a1bSJed Brown 
155*9371c9d4SSatish Balay     ii = 2;
156*9371c9d4SSatish Balay     jj = 4;
1579566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(preallocator, 1, &ii, 1, &jj, vv, INSERT_VALUES));
158c4762a1bSJed Brown 
159*9371c9d4SSatish Balay     ii = 4;
160*9371c9d4SSatish Balay     jj = 2;
1619566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(preallocator, 1, &ii, 1, &jj, vv, INSERT_VALUES));
162c4762a1bSJed Brown 
163*9371c9d4SSatish Balay     ii = 4;
164*9371c9d4SSatish Balay     jj = 4;
1659566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(preallocator, 1, &ii, 1, &jj, vv, INSERT_VALUES));
166c4762a1bSJed Brown 
167*9371c9d4SSatish Balay     ii = 9;
168*9371c9d4SSatish Balay     jj = 9;
1699566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(preallocator, 1, &ii, 1, &jj, vv, INSERT_VALUES));
170c4762a1bSJed Brown 
1719566063dSJacob Faibussowitsch     PetscCall(PetscFree(vv));
172c4762a1bSJed Brown   }
1739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
1749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   /*
177c4762a1bSJed Brown      Push non-zero pattern defined from preallocator into A.
178c4762a1bSJed Brown      Internally this will call MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE).
179c4762a1bSJed Brown      The arg fill = PETSC_TRUE will insert zeros in the matrix A automatically.
180c4762a1bSJed Brown   */
1819566063dSJacob Faibussowitsch   PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, A));
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   /*
184c4762a1bSJed Brown      We no longer require the preallocator object so destroy it.
185c4762a1bSJed Brown   */
1869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&preallocator));
187c4762a1bSJed Brown 
1889566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   {
191c4762a1bSJed Brown     PetscInt     ii, jj;
192c4762a1bSJed Brown     PetscScalar *vv;
193c4762a1bSJed Brown 
1949566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(bs * bs, &vv));
195*9371c9d4SSatish Balay     for (ii = 0; ii < bs * bs; ii++) { vv[ii] = (PetscReal)(ii + 1); }
196c4762a1bSJed Brown 
197*9371c9d4SSatish Balay     ii = 0;
198*9371c9d4SSatish Balay     jj = 9;
1999566063dSJacob Faibussowitsch     PetscCall(MatSetValue(A, ii, jj, 33.3, INSERT_VALUES));
200c4762a1bSJed Brown 
201*9371c9d4SSatish Balay     ii = 0;
202*9371c9d4SSatish Balay     jj = 0;
2039566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A, 1, &ii, 1, &jj, vv, INSERT_VALUES));
204c4762a1bSJed Brown 
205*9371c9d4SSatish Balay     ii = 2;
206*9371c9d4SSatish Balay     jj = 4;
2079566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A, 1, &ii, 1, &jj, vv, INSERT_VALUES));
208c4762a1bSJed Brown 
209*9371c9d4SSatish Balay     ii = 4;
210*9371c9d4SSatish Balay     jj = 2;
2119566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A, 1, &ii, 1, &jj, vv, INSERT_VALUES));
212c4762a1bSJed Brown 
213*9371c9d4SSatish Balay     ii = 4;
214*9371c9d4SSatish Balay     jj = 4;
2159566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A, 1, &ii, 1, &jj, vv, INSERT_VALUES));
216c4762a1bSJed Brown 
217*9371c9d4SSatish Balay     ii = 9;
218*9371c9d4SSatish Balay     jj = 9;
2199566063dSJacob Faibussowitsch     PetscCall(MatSetValuesBlocked(A, 1, &ii, 1, &jj, vv, INSERT_VALUES));
220c4762a1bSJed Brown 
2219566063dSJacob Faibussowitsch     PetscCall(PetscFree(vv));
222c4762a1bSJed Brown   }
2239566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2249566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
225c4762a1bSJed Brown 
2269566063dSJacob Faibussowitsch   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
227c4762a1bSJed Brown 
2289566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
229c4762a1bSJed Brown   PetscFunctionReturn(0);
230c4762a1bSJed Brown }
231c4762a1bSJed Brown 
232*9371c9d4SSatish Balay int main(int argc, char **args) {
233c4762a1bSJed Brown   PetscInt testid = 0;
234327415f7SBarry Smith   PetscFunctionBeginUser;
2359566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
2369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-test_id", &testid, NULL));
237c4762a1bSJed Brown   switch (testid) {
238*9371c9d4SSatish Balay   case 0: PetscCall(ex1_nonsquare_bs1()); break;
239*9371c9d4SSatish Balay   case 1: PetscCall(ex2_square_bsvariable()); break;
240*9371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Invalid value for -test_id. Must be {0,1}");
241c4762a1bSJed Brown   }
2429566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
243b122ec5aSJacob Faibussowitsch   return 0;
244c4762a1bSJed Brown }
245c4762a1bSJed Brown 
246c4762a1bSJed Brown /*TEST
247c4762a1bSJed Brown 
248c4762a1bSJed Brown    test:
249c4762a1bSJed Brown      suffix: t0_a_aij
250c4762a1bSJed Brown      nsize: 1
251c4762a1bSJed Brown      args: -test_id 0 -mat_type aij
252c4762a1bSJed Brown 
253c4762a1bSJed Brown    test:
254c4762a1bSJed Brown      suffix: t0_b_aij
255c4762a1bSJed Brown      nsize: 6
256c4762a1bSJed Brown      args: -test_id 0 -mat_type aij
257c4762a1bSJed Brown 
258c4762a1bSJed Brown    test:
259c4762a1bSJed Brown      suffix: t1_a_aij
260c4762a1bSJed Brown      nsize: 1
261c4762a1bSJed Brown      args: -test_id 1 -mat_type aij
262c4762a1bSJed Brown 
263c4762a1bSJed Brown    test:
264c4762a1bSJed Brown      suffix: t2_a_baij
265c4762a1bSJed Brown      nsize: 1
266c4762a1bSJed Brown      args: -test_id 1 -mat_type baij
267c4762a1bSJed Brown 
268c4762a1bSJed Brown    test:
269c4762a1bSJed Brown      suffix: t3_a_sbaij
270c4762a1bSJed Brown      nsize: 1
271c4762a1bSJed Brown      args: -test_id 1 -mat_type sbaij
272c4762a1bSJed Brown 
273c4762a1bSJed Brown    test:
274c4762a1bSJed Brown      suffix: t4_a_aij_bs3
275c4762a1bSJed Brown      nsize: 1
276c4762a1bSJed Brown      args: -test_id 1 -mat_type aij -block_size 3
277c4762a1bSJed Brown 
278c4762a1bSJed Brown    test:
279c4762a1bSJed Brown      suffix: t5_a_baij_bs3
280c4762a1bSJed Brown      nsize: 1
281c4762a1bSJed Brown      args: -test_id 1 -mat_type baij -block_size 3
282c4762a1bSJed Brown 
283c4762a1bSJed Brown    test:
284c4762a1bSJed Brown      suffix: t6_a_sbaij_bs3
285c4762a1bSJed Brown      nsize: 1
286c4762a1bSJed Brown      args: -test_id 1 -mat_type sbaij -block_size 3
287c4762a1bSJed Brown 
288c4762a1bSJed Brown    test:
289c4762a1bSJed Brown      suffix: t4_b_aij_bs3
290c4762a1bSJed Brown      nsize: 6
291c4762a1bSJed Brown      args: -test_id 1 -mat_type aij -block_size 3
292c4762a1bSJed Brown 
293c4762a1bSJed Brown    test:
294c4762a1bSJed Brown      suffix: t5_b_baij_bs3
295c4762a1bSJed Brown      nsize: 6
296c4762a1bSJed Brown      args: -test_id 1 -mat_type baij -block_size 3
297c4762a1bSJed Brown      filter: grep -v Mat_
298c4762a1bSJed Brown 
299c4762a1bSJed Brown    test:
300c4762a1bSJed Brown      suffix: t6_b_sbaij_bs3
301c4762a1bSJed Brown      nsize: 6
302c4762a1bSJed Brown      args: -test_id 1 -mat_type sbaij -block_size 3
303c4762a1bSJed Brown      filter: grep -v Mat_
304c4762a1bSJed Brown 
305c4762a1bSJed Brown TEST*/
306