xref: /petsc/src/mat/tests/ex5.c (revision 5f9962eec1b191700dd210fc035aa56572a01623)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatMult(), MatMultAdd(), MatMultTranspose().\n\
3*5f9962eeSHong Zhang Also MatMultTransposeAdd(), MatScale(), MatGetDiagonal(), MatDiagonalScale(), MatZeroEntries() and MatDuplicate().\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #include <petscmat.h>
6c4762a1bSJed Brown 
7d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
8d71ae5a4SJacob Faibussowitsch {
9c4762a1bSJed Brown   Mat         C;
10c4762a1bSJed Brown   Vec         s, u, w, x, y, z;
11c4762a1bSJed Brown   PetscInt    i, j, m = 8, n, rstart, rend, vstart, vend;
12c4762a1bSJed Brown   PetscScalar one = 1.0, negone = -1.0, v, alpha = 0.1;
13c4762a1bSJed Brown   PetscReal   norm, tol = PETSC_SQRT_MACHINE_EPSILON;
14c4762a1bSJed Brown   PetscBool   flg;
15c4762a1bSJed Brown 
16327415f7SBarry Smith   PetscFunctionBeginUser;
179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
189566063dSJacob Faibussowitsch   PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_COMMON));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
20c4762a1bSJed Brown   n = m;
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-rectA", &flg));
22c4762a1bSJed Brown   if (flg) n += 2;
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-rectB", &flg));
24c4762a1bSJed Brown   if (flg) n -= 2;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   /* ---------- Assemble matrix and vectors ----------- */
27c4762a1bSJed Brown 
289566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
299566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m, n));
309566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
319566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
329566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(C, &rstart, &rend));
339566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
349566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(x, PETSC_DECIDE, m));
359566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(x));
369566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &z));
379566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &w));
389566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &y));
399566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(y, PETSC_DECIDE, n));
409566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(y));
419566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(y, &u));
429566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(y, &s));
439566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(y, &vstart, &vend));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* Assembly */
46c4762a1bSJed Brown   for (i = rstart; i < rend; i++) {
47c4762a1bSJed Brown     v = 100 * (i + 1);
489566063dSJacob Faibussowitsch     PetscCall(VecSetValues(z, 1, &i, &v, INSERT_VALUES));
49c4762a1bSJed Brown     for (j = 0; j < n; j++) {
50c4762a1bSJed Brown       v = 10 * (i + 1) + j + 1;
519566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
52c4762a1bSJed Brown     }
53c4762a1bSJed Brown   }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /* Flush off proc Vec values and do more assembly */
569566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(z));
57c4762a1bSJed Brown   for (i = vstart; i < vend; i++) {
58c4762a1bSJed Brown     v = one * ((PetscReal)i);
599566063dSJacob Faibussowitsch     PetscCall(VecSetValues(y, 1, &i, &v, INSERT_VALUES));
60c4762a1bSJed Brown     v = 100.0 * i;
619566063dSJacob Faibussowitsch     PetscCall(VecSetValues(u, 1, &i, &v, INSERT_VALUES));
62c4762a1bSJed Brown   }
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Flush off proc Mat values and do more assembly */
659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FLUSH_ASSEMBLY));
66c4762a1bSJed Brown   for (i = rstart; i < rend; i++) {
67c4762a1bSJed Brown     for (j = 0; j < n; j++) {
68c4762a1bSJed Brown       v = 10 * (i + 1) + j + 1;
699566063dSJacob Faibussowitsch       PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
70c4762a1bSJed Brown     }
71c4762a1bSJed Brown   }
72c4762a1bSJed Brown   /* Try overlap Coomunication with the next stage XXXSetValues */
739566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(z));
74c4762a1bSJed Brown 
759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FLUSH_ASSEMBLY));
76c4762a1bSJed Brown   CHKMEMQ;
77c4762a1bSJed Brown   /* The Assembly for the second Stage */
789566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
799566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
809566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(y));
819566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(y));
829566063dSJacob Faibussowitsch   PetscCall(MatScale(C, alpha));
839566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(u));
849566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(u));
859566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "testing MatMult()\n"));
86c4762a1bSJed Brown   CHKMEMQ;
879566063dSJacob Faibussowitsch   PetscCall(MatMult(C, y, x));
88c4762a1bSJed Brown   CHKMEMQ;
899566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
909566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "testing MatMultAdd()\n"));
919566063dSJacob Faibussowitsch   PetscCall(MatMultAdd(C, y, z, w));
929566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x, one, z));
939566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x, negone, w));
949566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &norm));
9548a46eb9SPierre Jolivet   if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error difference = %g\n", (double)norm));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* ------- Test MatMultTranspose(), MatMultTransposeAdd() ------- */
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   for (i = rstart; i < rend; i++) {
100c4762a1bSJed Brown     v = one * ((PetscReal)i);
1019566063dSJacob Faibussowitsch     PetscCall(VecSetValues(x, 1, &i, &v, INSERT_VALUES));
102c4762a1bSJed Brown   }
1039566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(x));
1049566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(x));
1059566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "testing MatMultTranspose()\n"));
1069566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(C, x, y));
1079566063dSJacob Faibussowitsch   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
108c4762a1bSJed Brown 
1099566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "testing MatMultTransposeAdd()\n"));
1109566063dSJacob Faibussowitsch   PetscCall(MatMultTransposeAdd(C, x, u, s));
1119566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, one, u));
1129566063dSJacob Faibussowitsch   PetscCall(VecAXPY(y, negone, s));
1139566063dSJacob Faibussowitsch   PetscCall(VecNorm(y, NORM_2, &norm));
11448a46eb9SPierre Jolivet   if (norm > tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error difference = %g\n", (double)norm));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* -------------------- Test MatGetDiagonal() ------------------ */
117c4762a1bSJed Brown 
1189566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "testing MatGetDiagonal(), MatDiagonalScale()\n"));
1199566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
1209566063dSJacob Faibussowitsch   PetscCall(VecSet(x, one));
1219566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(C, x));
1229566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
123c4762a1bSJed Brown   for (i = vstart; i < vend; i++) {
124c4762a1bSJed Brown     v = one * ((PetscReal)(i + 1));
1259566063dSJacob Faibussowitsch     PetscCall(VecSetValues(y, 1, &i, &v, INSERT_VALUES));
126c4762a1bSJed Brown   }
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   /* -------------------- Test () MatDiagonalScale ------------------ */
1299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_diagonalscale", &flg));
130c4762a1bSJed Brown   if (flg) {
1319566063dSJacob Faibussowitsch     PetscCall(MatDiagonalScale(C, x, y));
1329566063dSJacob Faibussowitsch     PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
133c4762a1bSJed Brown   }
134*5f9962eeSHong Zhang   /* -------------------- Test () MatZeroEntries() and MatDuplicate() ------------------ */
135*5f9962eeSHong Zhang   PetscCall(PetscOptionsHasName(NULL, NULL, "-test_zeroentries", &flg));
136*5f9962eeSHong Zhang   if (flg) {
137*5f9962eeSHong Zhang     Mat D;
138*5f9962eeSHong Zhang     PetscCall(MatDuplicate(C, MAT_COPY_VALUES, &D));
139*5f9962eeSHong Zhang     PetscCall(MatZeroEntries(D));
140*5f9962eeSHong Zhang     PetscCall(MatView(D, PETSC_VIEWER_STDOUT_WORLD));
141*5f9962eeSHong Zhang     PetscCall(MatDestroy(&D));
142*5f9962eeSHong Zhang   }
143c4762a1bSJed Brown   /* Free data structures */
1449371c9d4SSatish Balay   PetscCall(VecDestroy(&u));
1459371c9d4SSatish Balay   PetscCall(VecDestroy(&s));
1469371c9d4SSatish Balay   PetscCall(VecDestroy(&w));
1479371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
1489371c9d4SSatish Balay   PetscCall(VecDestroy(&y));
1499371c9d4SSatish Balay   PetscCall(VecDestroy(&z));
1509566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
151c4762a1bSJed Brown 
1529566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
153b122ec5aSJacob Faibussowitsch   return 0;
154c4762a1bSJed Brown }
155c4762a1bSJed Brown 
156c4762a1bSJed Brown /*TEST
157c4762a1bSJed Brown 
158c4762a1bSJed Brown    test:
159c4762a1bSJed Brown       suffix: 11_A
160c4762a1bSJed Brown       args: -mat_type seqaij -rectA
161c4762a1bSJed Brown       filter: grep -v type
162c4762a1bSJed Brown 
163c4762a1bSJed Brown    test:
164c4762a1bSJed Brown       args: -mat_type seqdense -rectA
165c4762a1bSJed Brown       suffix: 12_A
166c4762a1bSJed Brown 
167c4762a1bSJed Brown    test:
168c4762a1bSJed Brown       args: -mat_type seqaij -rectB
169c4762a1bSJed Brown       suffix: 11_B
170c4762a1bSJed Brown       filter: grep -v type
171c4762a1bSJed Brown 
172c4762a1bSJed Brown    test:
173c4762a1bSJed Brown       args: -mat_type seqdense -rectB
174c4762a1bSJed Brown       suffix: 12_B
175c4762a1bSJed Brown 
176c4762a1bSJed Brown    test:
177c4762a1bSJed Brown       suffix: 21
178c4762a1bSJed Brown       args: -mat_type mpiaij
179c4762a1bSJed Brown       filter: grep -v type
180c4762a1bSJed Brown 
181c4762a1bSJed Brown    test:
182c4762a1bSJed Brown       suffix: 22
183c4762a1bSJed Brown       args: -mat_type mpidense
184c4762a1bSJed Brown 
185c4762a1bSJed Brown    test:
186c4762a1bSJed Brown       suffix: 23
187c4762a1bSJed Brown       nsize: 3
188c4762a1bSJed Brown       args: -mat_type mpiaij
189c4762a1bSJed Brown       filter: grep -v type
190c4762a1bSJed Brown 
191c4762a1bSJed Brown    test:
192c4762a1bSJed Brown       suffix: 24
193c4762a1bSJed Brown       nsize: 3
194c4762a1bSJed Brown       args: -mat_type mpidense
195c4762a1bSJed Brown 
196c4762a1bSJed Brown    test:
197c4762a1bSJed Brown       suffix: 2_aijcusparse_1
198c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda
199c4762a1bSJed Brown       filter: grep -v type
200c4762a1bSJed Brown       output_file: output/ex5_21.out
201c4762a1bSJed Brown       requires: cuda
202c4762a1bSJed Brown 
203c4762a1bSJed Brown    test:
204c4762a1bSJed Brown       nsize: 3
205c4762a1bSJed Brown       suffix: 2_aijcusparse_2
206c4762a1bSJed Brown       filter: grep -v type
207c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda
208bd46da1dSJunchao Zhang       args: -sf_type {{basic neighbor}}
209c4762a1bSJed Brown       output_file: output/ex5_23.out
210c4762a1bSJed Brown       requires: cuda
211c4762a1bSJed Brown 
212c4762a1bSJed Brown    test:
213c4762a1bSJed Brown       nsize: 3
214c4762a1bSJed Brown       suffix: 2_aijcusparse_3
215c4762a1bSJed Brown       filter: grep -v type
216c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda
217c20d7725SJed Brown       args: -sf_type {{basic neighbor}}
218c4762a1bSJed Brown       output_file: output/ex5_23.out
219dfd57a17SPierre Jolivet       requires: cuda defined(PETSC_HAVE_MPI_GPU_AWARE)
220c4762a1bSJed Brown 
221c4762a1bSJed Brown    test:
222c4762a1bSJed Brown       suffix: 31
223c4762a1bSJed Brown       args: -mat_type mpiaij -test_diagonalscale
224c4762a1bSJed Brown       filter: grep -v type
225c4762a1bSJed Brown 
226c4762a1bSJed Brown    test:
227c4762a1bSJed Brown       suffix: 32
228c4762a1bSJed Brown       args: -mat_type mpibaij -test_diagonalscale
229c4762a1bSJed Brown 
230c4762a1bSJed Brown    test:
231c4762a1bSJed Brown       suffix: 33
232c4762a1bSJed Brown       nsize: 3
233c4762a1bSJed Brown       args: -mat_type mpiaij -test_diagonalscale
234c4762a1bSJed Brown       filter: grep -v type
235c4762a1bSJed Brown 
236c4762a1bSJed Brown    test:
237c4762a1bSJed Brown       suffix: 34
238c4762a1bSJed Brown       nsize: 3
239c4762a1bSJed Brown       args: -mat_type mpibaij -test_diagonalscale
240c4762a1bSJed Brown 
241c4762a1bSJed Brown    test:
242c4762a1bSJed Brown       suffix: 3_aijcusparse_1
243c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda -test_diagonalscale
244c4762a1bSJed Brown       filter: grep -v type
245c4762a1bSJed Brown       output_file: output/ex5_31.out
246c4762a1bSJed Brown       requires: cuda
247c4762a1bSJed Brown 
248c4762a1bSJed Brown    test:
249c4762a1bSJed Brown       suffix: 3_aijcusparse_2
250c4762a1bSJed Brown       nsize: 3
251c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda -test_diagonalscale
252c4762a1bSJed Brown       filter: grep -v type
253c4762a1bSJed Brown       output_file: output/ex5_33.out
254c4762a1bSJed Brown       requires: cuda
255c4762a1bSJed Brown 
256c4762a1bSJed Brown    test:
25735990778SJunchao Zhang       suffix: 3_kokkos
25835990778SJunchao Zhang       nsize: 3
25935990778SJunchao Zhang       args: -mat_type mpiaijkokkos -vec_type kokkos -test_diagonalscale
26035990778SJunchao Zhang       filter: grep -v type
26135990778SJunchao Zhang       output_file: output/ex5_33.out
262dcfd994dSJunchao Zhang       requires: kokkos_kernels
26335990778SJunchao Zhang 
26435990778SJunchao Zhang    test:
265c4762a1bSJed Brown       suffix: aijcusparse_1
266c4762a1bSJed Brown       args: -mat_type seqaijcusparse -vec_type cuda -rectA
267c4762a1bSJed Brown       filter: grep -v type
268c4762a1bSJed Brown       output_file: output/ex5_11_A.out
269c4762a1bSJed Brown       requires: cuda
270c4762a1bSJed Brown 
271c4762a1bSJed Brown    test:
272c4762a1bSJed Brown       suffix: aijcusparse_2
273c4762a1bSJed Brown       args: -mat_type seqaijcusparse -vec_type cuda -rectB
274c4762a1bSJed Brown       filter: grep -v type
275c4762a1bSJed Brown       output_file: output/ex5_11_B.out
276c4762a1bSJed Brown       requires: cuda
277c4762a1bSJed Brown 
278c4762a1bSJed Brown    test:
279c4762a1bSJed Brown       suffix: sell_1
280*5f9962eeSHong Zhang       args: -mat_type sell -mat_sell_slice_height 8
281c4762a1bSJed Brown       output_file: output/ex5_41.out
282c4762a1bSJed Brown 
283c4762a1bSJed Brown    test:
284c4762a1bSJed Brown       suffix: sell_2
285c4762a1bSJed Brown       nsize: 3
286*5f9962eeSHong Zhang       args: -mat_type sell -mat_sell_slice_height 8
287c4762a1bSJed Brown       output_file: output/ex5_43.out
288c4762a1bSJed Brown 
289c4762a1bSJed Brown    test:
290c4762a1bSJed Brown       suffix: sell_3
291*5f9962eeSHong Zhang       args: -mat_type sell -test_diagonalscale -mat_sell_slice_height 8
292c4762a1bSJed Brown       output_file: output/ex5_51.out
293c4762a1bSJed Brown 
294c4762a1bSJed Brown    test:
295c4762a1bSJed Brown       suffix: sell_4
296c4762a1bSJed Brown       nsize: 3
297*5f9962eeSHong Zhang       args: -mat_type sell -test_diagonalscale -mat_sell_slice_height 8
298c4762a1bSJed Brown       output_file: output/ex5_53.out
299c4762a1bSJed Brown 
3002d1451d4SHong Zhang    test:
3012d1451d4SHong Zhang       suffix: sell_5
30290d2215bSHong Zhang       nsize: 3
303*5f9962eeSHong Zhang       args: -mat_type sellcuda -vec_type cuda -test_diagonalscale -test_zeroentries
30490d2215bSHong Zhang       output_file: output/ex5_55.out
3058711c661SHong Zhang       requires: cuda !complex
3062d1451d4SHong Zhang 
307*5f9962eeSHong Zhang    test:
308*5f9962eeSHong Zhang       suffix: sell_6
309*5f9962eeSHong Zhang       nsize: 3
310*5f9962eeSHong Zhang       args: -mat_type sellcuda -vec_type cuda -mat_sell_spmv_cuda_kernel {{1 2 3 4 5 6}}
311*5f9962eeSHong Zhang       output_file: output/ex5_56.out
312*5f9962eeSHong Zhang       requires: cuda !complex
313*5f9962eeSHong Zhang 
314*5f9962eeSHong Zhang    test:
315*5f9962eeSHong Zhang       suffix: sell_7
316*5f9962eeSHong Zhang       args: -m 32 -mat_type sellcuda -vec_type cuda -mat_sell_spmv_cuda_kernel {{0 7 9}} -mat_sell_spmv_cuda_blocky {{2 4 8 16 32}}
317*5f9962eeSHong Zhang       output_file: output/ex5_57.out
318*5f9962eeSHong Zhang       requires: cuda !complex !single
319c4762a1bSJed Brown TEST*/
320