xref: /petsc/src/mat/tests/ex230.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown static char help[] = "Example of using MatPreallocator\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*T
4c4762a1bSJed Brown    Concepts: Mat^matrix preallocation
5c4762a1bSJed Brown    Processors: n
6c4762a1bSJed Brown T*/
7c4762a1bSJed Brown 
8c4762a1bSJed Brown #include <petscmat.h>
9c4762a1bSJed Brown 
10c4762a1bSJed Brown PetscErrorCode ex1_nonsquare_bs1(void)
11c4762a1bSJed Brown {
12c4762a1bSJed Brown   Mat            A,preallocator;
13c4762a1bSJed Brown   PetscInt       M,N,m,n,bs;
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   /*
16c4762a1bSJed Brown      Create the Jacobian matrix
17c4762a1bSJed Brown   */
18362febeeSStefano Zampini   PetscFunctionBegin;
19c4762a1bSJed Brown   M = 10;
20c4762a1bSJed Brown   N = 8;
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATAIJ));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSize(A,1));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   /*
28c4762a1bSJed Brown      Get the sizes of the jacobian.
29c4762a1bSJed Brown   */
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(A,&m,&n));
31*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetBlockSize(A,&bs));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   /*
34c4762a1bSJed Brown      Create a preallocator matrix with sizes (local and global) matching the jacobian A.
35c4762a1bSJed Brown   */
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PetscObjectComm((PetscObject)A),&preallocator));
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(preallocator,MATPREALLOCATOR));
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(preallocator,m,n,M,N));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSize(preallocator,bs));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(preallocator));
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   /*
43c4762a1bSJed Brown      Insert non-zero pattern (e.g. perform a sweep over the grid).
44c4762a1bSJed Brown      You can use MatSetValues(), MatSetValuesBlocked() or MatSetValue().
45c4762a1bSJed Brown   */
46c4762a1bSJed Brown   {
47c4762a1bSJed Brown     PetscInt    ii,jj;
48c4762a1bSJed Brown     PetscScalar vv = 0.0;
49c4762a1bSJed Brown 
50c4762a1bSJed Brown     ii = 3; jj = 3;
51*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(preallocator,1,&ii,1,&jj,&vv,INSERT_VALUES));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown     ii = 7; jj = 4;
54*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(preallocator,1,&ii,1,&jj,&vv,INSERT_VALUES));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown     ii = 9; jj = 7;
57*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValue(preallocator,ii,jj,vv,INSERT_VALUES));
58c4762a1bSJed Brown   }
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY));
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /*
63c4762a1bSJed Brown      Push the non-zero pattern defined within preallocator into A.
64c4762a1bSJed Brown      Internally this will call MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE).
65c4762a1bSJed Brown      The arg fill = PETSC_TRUE will insert zeros in the matrix A automatically.
66c4762a1bSJed Brown   */
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /*
70c4762a1bSJed Brown      We no longer require the preallocator object so destroy it.
71c4762a1bSJed Brown   */
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&preallocator));
73c4762a1bSJed Brown 
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /*
77c4762a1bSJed Brown      Insert non-zero values into A.
78c4762a1bSJed Brown   */
79c4762a1bSJed Brown   {
80c4762a1bSJed Brown     PetscInt    ii,jj;
81c4762a1bSJed Brown     PetscScalar vv;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown     ii = 3; jj = 3; vv = 0.3;
84*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValue(A,ii,jj,vv,INSERT_VALUES));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown     ii = 7; jj = 4; vv = 3.3;
87*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&ii,1,&jj,&vv,INSERT_VALUES));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown     ii = 9; jj = 7; vv = 4.3;
90*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&ii,1,&jj,&vv,INSERT_VALUES));
91c4762a1bSJed Brown   }
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
94c4762a1bSJed Brown 
95*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
96c4762a1bSJed Brown 
97*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
98c4762a1bSJed Brown   PetscFunctionReturn(0);
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown PetscErrorCode ex2_square_bsvariable(void)
102c4762a1bSJed Brown {
103c4762a1bSJed Brown   Mat            A,preallocator;
104c4762a1bSJed Brown   PetscInt       M,N,m,n,bs = 1;
105c4762a1bSJed Brown 
106362febeeSStefano Zampini   PetscFunctionBegin;
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-block_size",&bs,NULL));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /*
110c4762a1bSJed Brown      Create the Jacobian matrix.
111c4762a1bSJed Brown   */
112c4762a1bSJed Brown   M = 10 * bs;
113c4762a1bSJed Brown   N = 10 * bs;
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATAIJ));
116*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSize(A,bs));
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /*
121c4762a1bSJed Brown      Get the sizes of the jacobian.
122c4762a1bSJed Brown   */
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(A,&m,&n));
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetBlockSize(A,&bs));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /*
127c4762a1bSJed Brown      Create a preallocator matrix with dimensions matching the jacobian A.
128c4762a1bSJed Brown   */
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PetscObjectComm((PetscObject)A),&preallocator));
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(preallocator,MATPREALLOCATOR));
131*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(preallocator,m,n,M,N));
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSize(preallocator,bs));
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(preallocator));
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /*
136c4762a1bSJed Brown      Insert non-zero pattern (e.g. perform a sweep over the grid).
137c4762a1bSJed Brown      You can use MatSetValues(), MatSetValuesBlocked() or MatSetValue().
138c4762a1bSJed Brown   */
139c4762a1bSJed Brown   {
140c4762a1bSJed Brown     PetscInt  ii,jj;
141c4762a1bSJed Brown     PetscScalar *vv;
142c4762a1bSJed Brown 
143*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(bs*bs,&vv));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown     ii = 0; jj = 9;
146*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValue(preallocator,ii,jj,vv[0],INSERT_VALUES));
147c4762a1bSJed Brown 
148c4762a1bSJed Brown     ii = 0; jj = 0;
149*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown     ii = 2; jj = 4;
152*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES));
153c4762a1bSJed Brown 
154c4762a1bSJed Brown     ii = 4; jj = 2;
155*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES));
156c4762a1bSJed Brown 
157c4762a1bSJed Brown     ii = 4; jj = 4;
158*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES));
159c4762a1bSJed Brown 
160c4762a1bSJed Brown     ii = 9; jj = 9;
161*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValuesBlocked(preallocator,1,&ii,1,&jj,vv,INSERT_VALUES));
162c4762a1bSJed Brown 
163*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(vv));
164c4762a1bSJed Brown   }
165*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(preallocator,MAT_FINAL_ASSEMBLY));
166*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(preallocator,MAT_FINAL_ASSEMBLY));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /*
169c4762a1bSJed Brown      Push non-zero pattern defined from preallocator into A.
170c4762a1bSJed Brown      Internally this will call MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE).
171c4762a1bSJed Brown      The arg fill = PETSC_TRUE will insert zeros in the matrix A automatically.
172c4762a1bSJed Brown   */
173*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPreallocatorPreallocate(preallocator,PETSC_TRUE,A));
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /*
176c4762a1bSJed Brown      We no longer require the preallocator object so destroy it.
177c4762a1bSJed Brown   */
178*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&preallocator));
179c4762a1bSJed Brown 
180*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   {
183c4762a1bSJed Brown     PetscInt    ii,jj;
184c4762a1bSJed Brown     PetscScalar *vv;
185c4762a1bSJed Brown 
186*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(bs*bs,&vv));
187c4762a1bSJed Brown     for (ii=0; ii<bs*bs; ii++) {
188c4762a1bSJed Brown       vv[ii] = (PetscReal)(ii+1);
189c4762a1bSJed Brown     }
190c4762a1bSJed Brown 
191c4762a1bSJed Brown     ii = 0; jj = 9;
192*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValue(A,ii,jj,33.3,INSERT_VALUES));
193c4762a1bSJed Brown 
194c4762a1bSJed Brown     ii = 0; jj = 0;
195*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES));
196c4762a1bSJed Brown 
197c4762a1bSJed Brown     ii = 2; jj = 4;
198*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES));
199c4762a1bSJed Brown 
200c4762a1bSJed Brown     ii = 4; jj = 2;
201*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown     ii = 4; jj = 4;
204*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES));
205c4762a1bSJed Brown 
206c4762a1bSJed Brown     ii = 9; jj = 9;
207*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValuesBlocked(A,1,&ii,1,&jj,vv,INSERT_VALUES));
208c4762a1bSJed Brown 
209*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(vv));
210c4762a1bSJed Brown   }
211*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
212*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
213c4762a1bSJed Brown 
214*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
215c4762a1bSJed Brown 
216*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
217c4762a1bSJed Brown   PetscFunctionReturn(0);
218c4762a1bSJed Brown }
219c4762a1bSJed Brown 
220c4762a1bSJed Brown int main(int argc, char **args)
221c4762a1bSJed Brown {
222c4762a1bSJed Brown   PetscErrorCode ierr;
223c4762a1bSJed Brown   PetscInt testid = 0;
224c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
225*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-test_id",&testid,NULL));
226c4762a1bSJed Brown   switch (testid) {
227c4762a1bSJed Brown     case 0:
228*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ex1_nonsquare_bs1());
229c4762a1bSJed Brown       break;
230c4762a1bSJed Brown     case 1:
231*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ex2_square_bsvariable());
232c4762a1bSJed Brown       break;
233c4762a1bSJed Brown     default:
234c4762a1bSJed Brown       SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Invalid value for -test_id. Must be {0,1}");
235c4762a1bSJed Brown   }
236c4762a1bSJed Brown   ierr = PetscFinalize();
237c4762a1bSJed Brown   return ierr;
238c4762a1bSJed Brown }
239c4762a1bSJed Brown 
240c4762a1bSJed Brown /*TEST
241c4762a1bSJed Brown 
242c4762a1bSJed Brown    test:
243c4762a1bSJed Brown      suffix: t0_a_aij
244c4762a1bSJed Brown      nsize: 1
245c4762a1bSJed Brown      args: -test_id 0 -mat_type aij
246c4762a1bSJed Brown 
247c4762a1bSJed Brown    test:
248c4762a1bSJed Brown      suffix: t0_b_aij
249c4762a1bSJed Brown      nsize: 6
250c4762a1bSJed Brown      args: -test_id 0 -mat_type aij
251c4762a1bSJed Brown 
252c4762a1bSJed Brown    test:
253c4762a1bSJed Brown      suffix: t1_a_aij
254c4762a1bSJed Brown      nsize: 1
255c4762a1bSJed Brown      args: -test_id 1 -mat_type aij
256c4762a1bSJed Brown 
257c4762a1bSJed Brown    test:
258c4762a1bSJed Brown      suffix: t2_a_baij
259c4762a1bSJed Brown      nsize: 1
260c4762a1bSJed Brown      args: -test_id 1 -mat_type baij
261c4762a1bSJed Brown 
262c4762a1bSJed Brown    test:
263c4762a1bSJed Brown      suffix: t3_a_sbaij
264c4762a1bSJed Brown      nsize: 1
265c4762a1bSJed Brown      args: -test_id 1 -mat_type sbaij
266c4762a1bSJed Brown 
267c4762a1bSJed Brown    test:
268c4762a1bSJed Brown      suffix: t4_a_aij_bs3
269c4762a1bSJed Brown      nsize: 1
270c4762a1bSJed Brown      args: -test_id 1 -mat_type aij -block_size 3
271c4762a1bSJed Brown 
272c4762a1bSJed Brown    test:
273c4762a1bSJed Brown      suffix: t5_a_baij_bs3
274c4762a1bSJed Brown      nsize: 1
275c4762a1bSJed Brown      args: -test_id 1 -mat_type baij -block_size 3
276c4762a1bSJed Brown 
277c4762a1bSJed Brown    test:
278c4762a1bSJed Brown      suffix: t6_a_sbaij_bs3
279c4762a1bSJed Brown      nsize: 1
280c4762a1bSJed Brown      args: -test_id 1 -mat_type sbaij -block_size 3
281c4762a1bSJed Brown 
282c4762a1bSJed Brown    test:
283c4762a1bSJed Brown      suffix: t4_b_aij_bs3
284c4762a1bSJed Brown      nsize: 6
285c4762a1bSJed Brown      args: -test_id 1 -mat_type aij -block_size 3
286c4762a1bSJed Brown 
287c4762a1bSJed Brown    test:
288c4762a1bSJed Brown      suffix: t5_b_baij_bs3
289c4762a1bSJed Brown      nsize: 6
290c4762a1bSJed Brown      args: -test_id 1 -mat_type baij -block_size 3
291c4762a1bSJed Brown      filter: grep -v Mat_
292c4762a1bSJed Brown 
293c4762a1bSJed Brown    test:
294c4762a1bSJed Brown      suffix: t6_b_sbaij_bs3
295c4762a1bSJed Brown      nsize: 6
296c4762a1bSJed Brown      args: -test_id 1 -mat_type sbaij -block_size 3
297c4762a1bSJed Brown      filter: grep -v Mat_
298c4762a1bSJed Brown 
299c4762a1bSJed Brown TEST*/
300