1c4762a1bSJed Brown static char help[] = "Example of using MatPreallocator\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown
ex1_nonsquare_bs1(void)5d71ae5a4SJacob Faibussowitsch PetscErrorCode ex1_nonsquare_bs1(void)
6d71ae5a4SJacob Faibussowitsch {
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
459371c9d4SSatish Balay ii = 3;
469371c9d4SSatish Balay jj = 3;
479566063dSJacob Faibussowitsch PetscCall(MatSetValues(preallocator, 1, &ii, 1, &jj, &vv, INSERT_VALUES));
48c4762a1bSJed Brown
499371c9d4SSatish Balay ii = 7;
509371c9d4SSatish Balay jj = 4;
519566063dSJacob Faibussowitsch PetscCall(MatSetValues(preallocator, 1, &ii, 1, &jj, &vv, INSERT_VALUES));
52c4762a1bSJed Brown
539371c9d4SSatish Balay ii = 9;
549371c9d4SSatish Balay jj = 7;
559566063dSJacob Faibussowitsch PetscCall(MatSetValue(preallocator, ii, jj, vv, INSERT_VALUES));
56c4762a1bSJed Brown }
579566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
589566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
59c4762a1bSJed Brown
60c4762a1bSJed Brown /*
61c4762a1bSJed Brown Push the non-zero pattern defined within preallocator into A.
62c4762a1bSJed Brown Internally this will call MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE).
63c4762a1bSJed Brown The arg fill = PETSC_TRUE will insert zeros in the matrix A automatically.
64c4762a1bSJed Brown */
659566063dSJacob Faibussowitsch PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, A));
66c4762a1bSJed Brown
67c4762a1bSJed Brown /*
68c4762a1bSJed Brown We no longer require the preallocator object so destroy it.
69c4762a1bSJed Brown */
709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&preallocator));
71c4762a1bSJed Brown
729566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
73c4762a1bSJed Brown
74c4762a1bSJed Brown /*
75c4762a1bSJed Brown Insert non-zero values into A.
76c4762a1bSJed Brown */
77c4762a1bSJed Brown {
78c4762a1bSJed Brown PetscInt ii, jj;
79c4762a1bSJed Brown PetscScalar vv;
80c4762a1bSJed Brown
819371c9d4SSatish Balay ii = 3;
829371c9d4SSatish Balay jj = 3;
839371c9d4SSatish Balay vv = 0.3;
849566063dSJacob Faibussowitsch PetscCall(MatSetValue(A, ii, jj, vv, INSERT_VALUES));
85c4762a1bSJed Brown
869371c9d4SSatish Balay ii = 7;
879371c9d4SSatish Balay jj = 4;
889371c9d4SSatish Balay vv = 3.3;
899566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &ii, 1, &jj, &vv, INSERT_VALUES));
90c4762a1bSJed Brown
919371c9d4SSatish Balay ii = 9;
929371c9d4SSatish Balay jj = 7;
939371c9d4SSatish Balay vv = 4.3;
949566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &ii, 1, &jj, &vv, INSERT_VALUES));
95c4762a1bSJed Brown }
969566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
979566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
98c4762a1bSJed Brown
999566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
100c4762a1bSJed Brown
1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
103c4762a1bSJed Brown }
104c4762a1bSJed Brown
ex2_square_bsvariable(void)105d71ae5a4SJacob Faibussowitsch PetscErrorCode ex2_square_bsvariable(void)
106d71ae5a4SJacob Faibussowitsch {
107c4762a1bSJed Brown Mat A, preallocator;
108c4762a1bSJed Brown PetscInt M, N, m, n, bs = 1;
109c4762a1bSJed Brown
110362febeeSStefano Zampini PetscFunctionBegin;
1119566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-block_size", &bs, NULL));
112c4762a1bSJed Brown
113c4762a1bSJed Brown /*
114c4762a1bSJed Brown Create the Jacobian matrix.
115c4762a1bSJed Brown */
116c4762a1bSJed Brown M = 10 * bs;
117c4762a1bSJed Brown N = 10 * bs;
1189566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1199566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATAIJ));
1209566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, M, N));
1219566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(A, bs));
1229566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
123c4762a1bSJed Brown
124c4762a1bSJed Brown /*
125c4762a1bSJed Brown Get the sizes of the jacobian.
126c4762a1bSJed Brown */
1279566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n));
1289566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &bs));
129c4762a1bSJed Brown
130c4762a1bSJed Brown /*
131c4762a1bSJed Brown Create a preallocator matrix with dimensions matching the jacobian A.
132c4762a1bSJed Brown */
1339566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &preallocator));
1349566063dSJacob Faibussowitsch PetscCall(MatSetType(preallocator, MATPREALLOCATOR));
1359566063dSJacob Faibussowitsch PetscCall(MatSetSizes(preallocator, m, n, M, N));
1369566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(preallocator, bs));
1379566063dSJacob Faibussowitsch PetscCall(MatSetUp(preallocator));
138c4762a1bSJed Brown
139c4762a1bSJed Brown /*
140c4762a1bSJed Brown Insert non-zero pattern (e.g. perform a sweep over the grid).
141c4762a1bSJed Brown You can use MatSetValues(), MatSetValuesBlocked() or MatSetValue().
142c4762a1bSJed Brown */
143c4762a1bSJed Brown {
144c4762a1bSJed Brown PetscInt ii, jj;
145c4762a1bSJed Brown PetscScalar *vv;
146c4762a1bSJed Brown
1479566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(bs * bs, &vv));
148c4762a1bSJed Brown
1499371c9d4SSatish Balay ii = 0;
1509371c9d4SSatish Balay jj = 9;
1519566063dSJacob Faibussowitsch PetscCall(MatSetValue(preallocator, ii, jj, vv[0], INSERT_VALUES));
152c4762a1bSJed Brown
1539371c9d4SSatish Balay ii = 0;
1549371c9d4SSatish Balay jj = 0;
1559566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(preallocator, 1, &ii, 1, &jj, vv, INSERT_VALUES));
156c4762a1bSJed Brown
1579371c9d4SSatish Balay ii = 2;
1589371c9d4SSatish Balay jj = 4;
1599566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(preallocator, 1, &ii, 1, &jj, vv, INSERT_VALUES));
160c4762a1bSJed Brown
1619371c9d4SSatish Balay ii = 4;
1629371c9d4SSatish Balay jj = 2;
1639566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(preallocator, 1, &ii, 1, &jj, vv, INSERT_VALUES));
164c4762a1bSJed Brown
1659371c9d4SSatish Balay ii = 4;
1669371c9d4SSatish Balay jj = 4;
1679566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(preallocator, 1, &ii, 1, &jj, vv, INSERT_VALUES));
168c4762a1bSJed Brown
1699371c9d4SSatish Balay ii = 9;
1709371c9d4SSatish Balay jj = 9;
1719566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(preallocator, 1, &ii, 1, &jj, vv, INSERT_VALUES));
172c4762a1bSJed Brown
1739566063dSJacob Faibussowitsch PetscCall(PetscFree(vv));
174c4762a1bSJed Brown }
1759566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY));
1769566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY));
177c4762a1bSJed Brown
178c4762a1bSJed Brown /*
179c4762a1bSJed Brown Push non-zero pattern defined from preallocator into A.
180c4762a1bSJed Brown Internally this will call MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE).
181c4762a1bSJed Brown The arg fill = PETSC_TRUE will insert zeros in the matrix A automatically.
182c4762a1bSJed Brown */
1839566063dSJacob Faibussowitsch PetscCall(MatPreallocatorPreallocate(preallocator, PETSC_TRUE, A));
184c4762a1bSJed Brown
185c4762a1bSJed Brown /*
186c4762a1bSJed Brown We no longer require the preallocator object so destroy it.
187c4762a1bSJed Brown */
1889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&preallocator));
189c4762a1bSJed Brown
1909566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
191c4762a1bSJed Brown
192c4762a1bSJed Brown {
193c4762a1bSJed Brown PetscInt ii, jj;
194c4762a1bSJed Brown PetscScalar *vv;
195c4762a1bSJed Brown
1969566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(bs * bs, &vv));
197ad540459SPierre Jolivet for (ii = 0; ii < bs * bs; ii++) vv[ii] = (PetscReal)(ii + 1);
198c4762a1bSJed Brown
1999371c9d4SSatish Balay ii = 0;
2009371c9d4SSatish Balay jj = 9;
2019566063dSJacob Faibussowitsch PetscCall(MatSetValue(A, ii, jj, 33.3, INSERT_VALUES));
202c4762a1bSJed Brown
2039371c9d4SSatish Balay ii = 0;
2049371c9d4SSatish Balay jj = 0;
2059566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(A, 1, &ii, 1, &jj, vv, INSERT_VALUES));
206c4762a1bSJed Brown
2079371c9d4SSatish Balay ii = 2;
2089371c9d4SSatish Balay jj = 4;
2099566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(A, 1, &ii, 1, &jj, vv, INSERT_VALUES));
210c4762a1bSJed Brown
2119371c9d4SSatish Balay ii = 4;
2129371c9d4SSatish Balay jj = 2;
2139566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(A, 1, &ii, 1, &jj, vv, INSERT_VALUES));
214c4762a1bSJed Brown
2159371c9d4SSatish Balay ii = 4;
2169371c9d4SSatish Balay jj = 4;
2179566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(A, 1, &ii, 1, &jj, vv, INSERT_VALUES));
218c4762a1bSJed Brown
2199371c9d4SSatish Balay ii = 9;
2209371c9d4SSatish Balay jj = 9;
2219566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlocked(A, 1, &ii, 1, &jj, vv, INSERT_VALUES));
222c4762a1bSJed Brown
2239566063dSJacob Faibussowitsch PetscCall(PetscFree(vv));
224c4762a1bSJed Brown }
2259566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2269566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
227c4762a1bSJed Brown
2289566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
229c4762a1bSJed Brown
2309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
2313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
232c4762a1bSJed Brown }
233c4762a1bSJed Brown
main(int argc,char ** args)234d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
235d71ae5a4SJacob Faibussowitsch {
236c4762a1bSJed Brown PetscInt testid = 0;
2374d86920dSPierre Jolivet
238327415f7SBarry Smith PetscFunctionBeginUser;
239*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
2409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-test_id", &testid, NULL));
241c4762a1bSJed Brown switch (testid) {
242d71ae5a4SJacob Faibussowitsch case 0:
243d71ae5a4SJacob Faibussowitsch PetscCall(ex1_nonsquare_bs1());
244d71ae5a4SJacob Faibussowitsch break;
245d71ae5a4SJacob Faibussowitsch case 1:
246d71ae5a4SJacob Faibussowitsch PetscCall(ex2_square_bsvariable());
247d71ae5a4SJacob Faibussowitsch break;
248d71ae5a4SJacob Faibussowitsch default:
249d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Invalid value for -test_id. Must be {0,1}");
250c4762a1bSJed Brown }
2519566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
252b122ec5aSJacob Faibussowitsch return 0;
253c4762a1bSJed Brown }
254c4762a1bSJed Brown
255c4762a1bSJed Brown /*TEST
256c4762a1bSJed Brown
257c4762a1bSJed Brown test:
258c4762a1bSJed Brown suffix: t0_a_aij
259c4762a1bSJed Brown nsize: 1
260c4762a1bSJed Brown args: -test_id 0 -mat_type aij
261c4762a1bSJed Brown
262c4762a1bSJed Brown test:
263c4762a1bSJed Brown suffix: t0_b_aij
264c4762a1bSJed Brown nsize: 6
265c4762a1bSJed Brown args: -test_id 0 -mat_type aij
266c4762a1bSJed Brown
267c4762a1bSJed Brown test:
268c4762a1bSJed Brown suffix: t1_a_aij
269c4762a1bSJed Brown nsize: 1
270c4762a1bSJed Brown args: -test_id 1 -mat_type aij
271c4762a1bSJed Brown
272c4762a1bSJed Brown test:
273c4762a1bSJed Brown suffix: t2_a_baij
274c4762a1bSJed Brown nsize: 1
275c4762a1bSJed Brown args: -test_id 1 -mat_type baij
276c4762a1bSJed Brown
277c4762a1bSJed Brown test:
278c4762a1bSJed Brown suffix: t3_a_sbaij
279c4762a1bSJed Brown nsize: 1
280c4762a1bSJed Brown args: -test_id 1 -mat_type sbaij
281c4762a1bSJed Brown
282c4762a1bSJed Brown test:
283c4762a1bSJed Brown suffix: t4_a_aij_bs3
284c4762a1bSJed Brown nsize: 1
285c4762a1bSJed Brown args: -test_id 1 -mat_type aij -block_size 3
286c4762a1bSJed Brown
287c4762a1bSJed Brown test:
288c4762a1bSJed Brown suffix: t5_a_baij_bs3
289c4762a1bSJed Brown nsize: 1
290c4762a1bSJed Brown args: -test_id 1 -mat_type baij -block_size 3
291c4762a1bSJed Brown
292c4762a1bSJed Brown test:
293c4762a1bSJed Brown suffix: t6_a_sbaij_bs3
294c4762a1bSJed Brown nsize: 1
295c4762a1bSJed Brown args: -test_id 1 -mat_type sbaij -block_size 3
296c4762a1bSJed Brown
297c4762a1bSJed Brown test:
298c4762a1bSJed Brown suffix: t4_b_aij_bs3
299c4762a1bSJed Brown nsize: 6
300c4762a1bSJed Brown args: -test_id 1 -mat_type aij -block_size 3
301c4762a1bSJed Brown
302c4762a1bSJed Brown test:
303c4762a1bSJed Brown suffix: t5_b_baij_bs3
304c4762a1bSJed Brown nsize: 6
305c4762a1bSJed Brown args: -test_id 1 -mat_type baij -block_size 3
306c4762a1bSJed Brown
307c4762a1bSJed Brown test:
308c4762a1bSJed Brown suffix: t6_b_sbaij_bs3
309c4762a1bSJed Brown nsize: 6
310c4762a1bSJed Brown args: -test_id 1 -mat_type sbaij -block_size 3
311c4762a1bSJed Brown
312c4762a1bSJed Brown TEST*/
313