xref: /petsc/src/mat/tests/ex5.c (revision 3599077843b2bc2bea7ec9ce431e983cccba183a)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatMult(), MatMultAdd(), MatMultTranspose().\n\
3c4762a1bSJed Brown Also MatMultTransposeAdd(), MatScale(), MatGetDiagonal(), and MatDiagonalScale().\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #include <petscmat.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown int main(int argc,char **args)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   Mat            C;
10c4762a1bSJed Brown   Vec            s,u,w,x,y,z;
11c4762a1bSJed Brown   PetscErrorCode ierr;
12c4762a1bSJed Brown   PetscInt       i,j,m = 8,n,rstart,rend,vstart,vend;
13c4762a1bSJed Brown   PetscScalar    one = 1.0,negone = -1.0,v,alpha=0.1;
14c4762a1bSJed Brown   PetscReal      norm, tol = PETSC_SQRT_MACHINE_EPSILON;
15c4762a1bSJed Brown   PetscBool      flg;
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
18c4762a1bSJed Brown   ierr = PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON);CHKERRQ(ierr);
19c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
20c4762a1bSJed Brown   n    = m;
21c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-rectA",&flg);CHKERRQ(ierr);
22c4762a1bSJed Brown   if (flg) n += 2;
23c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-rectB",&flg);CHKERRQ(ierr);
24c4762a1bSJed Brown   if (flg) n -= 2;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   /* ---------- Assemble matrix and vectors ----------- */
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
29c4762a1bSJed Brown   ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr = MatSetUp(C);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr);
33c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
34c4762a1bSJed Brown   ierr = VecSetSizes(x,PETSC_DECIDE,m);CHKERRQ(ierr);
35c4762a1bSJed Brown   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
36c4762a1bSJed Brown   ierr = VecDuplicate(x,&z);CHKERRQ(ierr);
37c4762a1bSJed Brown   ierr = VecDuplicate(x,&w);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&y);CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = VecSetSizes(y,PETSC_DECIDE,n);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = VecSetFromOptions(y);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = VecDuplicate(y,&u);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = VecDuplicate(y,&s);CHKERRQ(ierr);
43c4762a1bSJed Brown   ierr = VecGetOwnershipRange(y,&vstart,&vend);CHKERRQ(ierr);
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Assembly */
46c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
47c4762a1bSJed Brown     v    = 100*(i+1);
48c4762a1bSJed Brown     ierr = VecSetValues(z,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr);
49c4762a1bSJed Brown     for (j=0; j<n; j++) {
50c4762a1bSJed Brown       v    = 10*(i+1)+j+1;
51c4762a1bSJed Brown       ierr = MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
52c4762a1bSJed Brown     }
53c4762a1bSJed Brown   }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /* Flush off proc Vec values and do more assembly */
56c4762a1bSJed Brown   ierr = VecAssemblyBegin(z);CHKERRQ(ierr);
57c4762a1bSJed Brown   for (i=vstart; i<vend; i++) {
58c4762a1bSJed Brown     v    = one*((PetscReal)i);
59c4762a1bSJed Brown     ierr = VecSetValues(y,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr);
60c4762a1bSJed Brown     v    = 100.0*i;
61c4762a1bSJed Brown     ierr = VecSetValues(u,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr);
62c4762a1bSJed Brown   }
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Flush off proc Mat values and do more assembly */
65c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
66c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
67c4762a1bSJed Brown     for (j=0; j<n; j++) {
68c4762a1bSJed Brown       v    = 10*(i+1)+j+1;
69c4762a1bSJed Brown       ierr = MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
70c4762a1bSJed Brown     }
71c4762a1bSJed Brown   }
72c4762a1bSJed Brown   /* Try overlap Coomunication with the next stage XXXSetValues */
73c4762a1bSJed Brown   ierr = VecAssemblyEnd(z);CHKERRQ(ierr);
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
76c4762a1bSJed Brown   CHKMEMQ;
77c4762a1bSJed Brown   /* The Assembly for the second Stage */
78c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
79c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
80c4762a1bSJed Brown   ierr = VecAssemblyBegin(y);CHKERRQ(ierr);
81c4762a1bSJed Brown   ierr = VecAssemblyEnd(y);CHKERRQ(ierr);
82c4762a1bSJed Brown   ierr = MatScale(C,alpha);CHKERRQ(ierr);
83c4762a1bSJed Brown   ierr = VecAssemblyBegin(u);CHKERRQ(ierr);
84c4762a1bSJed Brown   ierr = VecAssemblyEnd(u);CHKERRQ(ierr);
85c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"testing MatMult()\n");CHKERRQ(ierr);
86c4762a1bSJed Brown   CHKMEMQ;
87c4762a1bSJed Brown   ierr = MatMult(C,y,x);CHKERRQ(ierr);
88c4762a1bSJed Brown   CHKMEMQ;
89c4762a1bSJed Brown   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
90c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"testing MatMultAdd()\n");CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = MatMultAdd(C,y,z,w);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = VecAXPY(x,one,z);CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr = VecAXPY(x,negone,w);CHKERRQ(ierr);
94c4762a1bSJed Brown   ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr);
95c4762a1bSJed Brown   if (norm > tol) {
96c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm of error difference = %g\n",(double)norm);CHKERRQ(ierr);
97c4762a1bSJed Brown   }
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* ------- Test MatMultTranspose(), MatMultTransposeAdd() ------- */
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
102c4762a1bSJed Brown     v    = one*((PetscReal)i);
103c4762a1bSJed Brown     ierr = VecSetValues(x,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr);
104c4762a1bSJed Brown   }
105c4762a1bSJed Brown   ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
106c4762a1bSJed Brown   ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
107c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"testing MatMultTranspose()\n");CHKERRQ(ierr);
108c4762a1bSJed Brown   ierr = MatMultTranspose(C,x,y);CHKERRQ(ierr);
109c4762a1bSJed Brown   ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"testing MatMultTransposeAdd()\n");CHKERRQ(ierr);
112c4762a1bSJed Brown   ierr = MatMultTransposeAdd(C,x,u,s);CHKERRQ(ierr);
113c4762a1bSJed Brown   ierr = VecAXPY(y,one,u);CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr = VecAXPY(y,negone,s);CHKERRQ(ierr);
115c4762a1bSJed Brown   ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
116c4762a1bSJed Brown   if (norm > tol) {
117c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm of error difference = %g\n",(double)norm);CHKERRQ(ierr);
118c4762a1bSJed Brown   }
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /* -------------------- Test MatGetDiagonal() ------------------ */
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"testing MatGetDiagonal(), MatDiagonalScale()\n");CHKERRQ(ierr);
123c4762a1bSJed Brown   ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
124c4762a1bSJed Brown   ierr = VecSet(x,one);CHKERRQ(ierr);
125c4762a1bSJed Brown   ierr = MatGetDiagonal(C,x);CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
127c4762a1bSJed Brown   for (i=vstart; i<vend; i++) {
128c4762a1bSJed Brown     v    = one*((PetscReal)(i+1));
129c4762a1bSJed Brown     ierr = VecSetValues(y,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr);
130c4762a1bSJed Brown   }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* -------------------- Test () MatDiagonalScale ------------------ */
133c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-test_diagonalscale",&flg);CHKERRQ(ierr);
134c4762a1bSJed Brown   if (flg) {
135c4762a1bSJed Brown     ierr = MatDiagonalScale(C,x,y);CHKERRQ(ierr);
136c4762a1bSJed Brown     ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
137c4762a1bSJed Brown   }
138c4762a1bSJed Brown   /* Free data structures */
139c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr); ierr = VecDestroy(&s);CHKERRQ(ierr);
140c4762a1bSJed Brown   ierr = VecDestroy(&w);CHKERRQ(ierr); ierr = VecDestroy(&x);CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr = VecDestroy(&y);CHKERRQ(ierr); ierr = VecDestroy(&z);CHKERRQ(ierr);
142c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   ierr = PetscFinalize();
145c4762a1bSJed Brown   return ierr;
146c4762a1bSJed Brown }
147c4762a1bSJed Brown 
148c4762a1bSJed Brown 
149c4762a1bSJed Brown /*TEST
150c4762a1bSJed Brown 
151c4762a1bSJed Brown 
152c4762a1bSJed Brown    test:
153c4762a1bSJed Brown       suffix: 11_A
154c4762a1bSJed Brown       args: -mat_type seqaij -rectA
155c4762a1bSJed Brown       filter: grep -v type
156c4762a1bSJed Brown 
157c4762a1bSJed Brown    test:
158c4762a1bSJed Brown       args: -mat_type seqdense -rectA
159c4762a1bSJed Brown       suffix: 12_A
160c4762a1bSJed Brown 
161c4762a1bSJed Brown    test:
162c4762a1bSJed Brown       args: -mat_type seqaij -rectB
163c4762a1bSJed Brown       suffix: 11_B
164c4762a1bSJed Brown       filter: grep -v type
165c4762a1bSJed Brown 
166c4762a1bSJed Brown    test:
167c4762a1bSJed Brown       args: -mat_type seqdense -rectB
168c4762a1bSJed Brown       suffix: 12_B
169c4762a1bSJed Brown 
170c4762a1bSJed Brown    test:
171c4762a1bSJed Brown       suffix: 21
172c4762a1bSJed Brown       args: -mat_type mpiaij
173c4762a1bSJed Brown       filter: grep -v type
174c4762a1bSJed Brown 
175c4762a1bSJed Brown    test:
176c4762a1bSJed Brown       suffix: 22
177c4762a1bSJed Brown       args: -mat_type mpidense
178c4762a1bSJed Brown 
179c4762a1bSJed Brown    test:
180c4762a1bSJed Brown       suffix: 23
181c4762a1bSJed Brown       nsize: 3
182c4762a1bSJed Brown       args: -mat_type mpiaij
183c4762a1bSJed Brown       filter: grep -v type
184c4762a1bSJed Brown 
185c4762a1bSJed Brown    test:
186c4762a1bSJed Brown       suffix: 24
187c4762a1bSJed Brown       nsize: 3
188c4762a1bSJed Brown       args: -mat_type mpidense
189c4762a1bSJed Brown 
190c4762a1bSJed Brown    test:
191c4762a1bSJed Brown       suffix: 2_aijcusparse_1
192c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda
193c4762a1bSJed Brown       filter: grep -v type
194c4762a1bSJed Brown       output_file: output/ex5_21.out
195c4762a1bSJed Brown       requires: cuda
196c4762a1bSJed Brown 
197c4762a1bSJed Brown 
198c4762a1bSJed Brown    test:
199c4762a1bSJed Brown       nsize: 3
200c4762a1bSJed Brown       suffix: 2_aijcusparse_2
201c4762a1bSJed Brown       filter: grep -v type
202c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda
203c20d7725SJed Brown       args: -sf_type {{basic neighbor}} -vecscatter_packongpu {{0 1}}
204c4762a1bSJed Brown       output_file: output/ex5_23.out
205c4762a1bSJed Brown       requires: cuda
206c4762a1bSJed Brown 
207*35990778SJunchao Zhang 
208c4762a1bSJed Brown    test:
209c4762a1bSJed Brown       nsize: 3
210c4762a1bSJed Brown       suffix: 2_aijcusparse_3
211c4762a1bSJed Brown       filter: grep -v type
212c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda
213c20d7725SJed Brown       args: -sf_type {{basic neighbor}}
214c4762a1bSJed Brown       output_file: output/ex5_23.out
215c4762a1bSJed Brown       requires: cuda define(PETSC_HAVE_MPI_GPU_AWARE)
216c4762a1bSJed Brown 
217c4762a1bSJed Brown    test:
218c4762a1bSJed Brown       suffix: 31
219c4762a1bSJed Brown       args: -mat_type mpiaij -test_diagonalscale
220c4762a1bSJed Brown       filter: grep -v type
221c4762a1bSJed Brown 
222c4762a1bSJed Brown    test:
223c4762a1bSJed Brown       suffix: 32
224c4762a1bSJed Brown       args: -mat_type mpibaij -test_diagonalscale
225c4762a1bSJed Brown       filter: grep -v Mat_
226c4762a1bSJed Brown 
227c4762a1bSJed Brown    test:
228c4762a1bSJed Brown       suffix: 33
229c4762a1bSJed Brown       nsize: 3
230c4762a1bSJed Brown       args: -mat_type mpiaij -test_diagonalscale
231c4762a1bSJed Brown       filter: grep -v type
232c4762a1bSJed Brown 
233c4762a1bSJed Brown    test:
234c4762a1bSJed Brown       suffix: 34
235c4762a1bSJed Brown       nsize: 3
236c4762a1bSJed Brown       args: -mat_type mpibaij -test_diagonalscale
237c4762a1bSJed Brown       filter: grep -v Mat_
238c4762a1bSJed Brown 
239c4762a1bSJed Brown    test:
240c4762a1bSJed Brown       suffix: 3_aijcusparse_1
241c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda -test_diagonalscale
242c4762a1bSJed Brown       filter: grep -v type
243c4762a1bSJed Brown       output_file: output/ex5_31.out
244c4762a1bSJed Brown       requires: cuda
245c4762a1bSJed Brown 
246c4762a1bSJed Brown    test:
247c4762a1bSJed Brown       suffix: 3_aijcusparse_2
248c4762a1bSJed Brown       nsize: 3
249c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda -test_diagonalscale
250c4762a1bSJed Brown       filter: grep -v type
251c4762a1bSJed Brown       output_file: output/ex5_33.out
252c4762a1bSJed Brown       requires: cuda
253c4762a1bSJed Brown 
254c4762a1bSJed Brown    test:
255*35990778SJunchao Zhang       suffix: 3_kokkos
256*35990778SJunchao Zhang       nsize: 3
257*35990778SJunchao Zhang       args: -mat_type mpiaijkokkos -vec_type kokkos -test_diagonalscale
258*35990778SJunchao Zhang       filter: grep -v type
259*35990778SJunchao Zhang       output_file: output/ex5_33.out
260*35990778SJunchao Zhang       requires: kokkos_kernels
261*35990778SJunchao Zhang 
262*35990778SJunchao Zhang    test:
263c4762a1bSJed Brown       suffix: aijcusparse_1
264c4762a1bSJed Brown       args: -mat_type seqaijcusparse -vec_type cuda -rectA
265c4762a1bSJed Brown       filter: grep -v type
266c4762a1bSJed Brown       output_file: output/ex5_11_A.out
267c4762a1bSJed Brown       requires: cuda
268c4762a1bSJed Brown 
269c4762a1bSJed Brown    test:
270c4762a1bSJed Brown       suffix: aijcusparse_2
271c4762a1bSJed Brown       args: -mat_type seqaijcusparse -vec_type cuda -rectB
272c4762a1bSJed Brown       filter: grep -v type
273c4762a1bSJed Brown       output_file: output/ex5_11_B.out
274c4762a1bSJed Brown       requires: cuda
275c4762a1bSJed Brown 
276c4762a1bSJed Brown    test:
277c4762a1bSJed Brown       suffix: sell_1
278c4762a1bSJed Brown       args: -mat_type sell
279c4762a1bSJed Brown       output_file: output/ex5_41.out
280c4762a1bSJed Brown 
281c4762a1bSJed Brown    test:
282c4762a1bSJed Brown       suffix: sell_2
283c4762a1bSJed Brown       nsize: 3
284c4762a1bSJed Brown       args: -mat_type sell
285c4762a1bSJed Brown       output_file: output/ex5_43.out
286c4762a1bSJed Brown 
287c4762a1bSJed Brown    test:
288c4762a1bSJed Brown       suffix: sell_3
289c4762a1bSJed Brown       args: -mat_type sell -test_diagonalscale
290c4762a1bSJed Brown       output_file: output/ex5_51.out
291c4762a1bSJed Brown 
292c4762a1bSJed Brown    test:
293c4762a1bSJed Brown       suffix: sell_4
294c4762a1bSJed Brown       nsize: 3
295c4762a1bSJed Brown       args: -mat_type sell -test_diagonalscale
296c4762a1bSJed Brown       output_file: output/ex5_53.out
297c4762a1bSJed Brown 
298c4762a1bSJed Brown TEST*/
299