xref: /petsc/src/mat/tests/ex5.c (revision 8711c6618dd802f17a929eb273bae3d4da533dd7)
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 
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   }
134c4762a1bSJed Brown   /* Free data structures */
1359371c9d4SSatish Balay   PetscCall(VecDestroy(&u));
1369371c9d4SSatish Balay   PetscCall(VecDestroy(&s));
1379371c9d4SSatish Balay   PetscCall(VecDestroy(&w));
1389371c9d4SSatish Balay   PetscCall(VecDestroy(&x));
1399371c9d4SSatish Balay   PetscCall(VecDestroy(&y));
1409371c9d4SSatish Balay   PetscCall(VecDestroy(&z));
1419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
142c4762a1bSJed Brown 
1439566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
144b122ec5aSJacob Faibussowitsch   return 0;
145c4762a1bSJed Brown }
146c4762a1bSJed Brown 
147c4762a1bSJed Brown /*TEST
148c4762a1bSJed Brown 
149c4762a1bSJed Brown    test:
150c4762a1bSJed Brown       suffix: 11_A
151c4762a1bSJed Brown       args: -mat_type seqaij -rectA
152c4762a1bSJed Brown       filter: grep -v type
153c4762a1bSJed Brown 
154c4762a1bSJed Brown    test:
155c4762a1bSJed Brown       args: -mat_type seqdense -rectA
156c4762a1bSJed Brown       suffix: 12_A
157c4762a1bSJed Brown 
158c4762a1bSJed Brown    test:
159c4762a1bSJed Brown       args: -mat_type seqaij -rectB
160c4762a1bSJed Brown       suffix: 11_B
161c4762a1bSJed Brown       filter: grep -v type
162c4762a1bSJed Brown 
163c4762a1bSJed Brown    test:
164c4762a1bSJed Brown       args: -mat_type seqdense -rectB
165c4762a1bSJed Brown       suffix: 12_B
166c4762a1bSJed Brown 
167c4762a1bSJed Brown    test:
168c4762a1bSJed Brown       suffix: 21
169c4762a1bSJed Brown       args: -mat_type mpiaij
170c4762a1bSJed Brown       filter: grep -v type
171c4762a1bSJed Brown 
172c4762a1bSJed Brown    test:
173c4762a1bSJed Brown       suffix: 22
174c4762a1bSJed Brown       args: -mat_type mpidense
175c4762a1bSJed Brown 
176c4762a1bSJed Brown    test:
177c4762a1bSJed Brown       suffix: 23
178c4762a1bSJed Brown       nsize: 3
179c4762a1bSJed Brown       args: -mat_type mpiaij
180c4762a1bSJed Brown       filter: grep -v type
181c4762a1bSJed Brown 
182c4762a1bSJed Brown    test:
183c4762a1bSJed Brown       suffix: 24
184c4762a1bSJed Brown       nsize: 3
185c4762a1bSJed Brown       args: -mat_type mpidense
186c4762a1bSJed Brown 
187c4762a1bSJed Brown    test:
188c4762a1bSJed Brown       suffix: 2_aijcusparse_1
189c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda
190c4762a1bSJed Brown       filter: grep -v type
191c4762a1bSJed Brown       output_file: output/ex5_21.out
192c4762a1bSJed Brown       requires: cuda
193c4762a1bSJed Brown 
194c4762a1bSJed Brown    test:
195c4762a1bSJed Brown       nsize: 3
196c4762a1bSJed Brown       suffix: 2_aijcusparse_2
197c4762a1bSJed Brown       filter: grep -v type
198c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda
199bd46da1dSJunchao Zhang       args: -sf_type {{basic neighbor}}
200c4762a1bSJed Brown       output_file: output/ex5_23.out
201c4762a1bSJed Brown       requires: cuda
202c4762a1bSJed Brown 
203c4762a1bSJed Brown    test:
204c4762a1bSJed Brown       nsize: 3
205c4762a1bSJed Brown       suffix: 2_aijcusparse_3
206c4762a1bSJed Brown       filter: grep -v type
207c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda
208c20d7725SJed Brown       args: -sf_type {{basic neighbor}}
209c4762a1bSJed Brown       output_file: output/ex5_23.out
210dfd57a17SPierre Jolivet       requires: cuda defined(PETSC_HAVE_MPI_GPU_AWARE)
211c4762a1bSJed Brown 
212c4762a1bSJed Brown    test:
213c4762a1bSJed Brown       suffix: 31
214c4762a1bSJed Brown       args: -mat_type mpiaij -test_diagonalscale
215c4762a1bSJed Brown       filter: grep -v type
216c4762a1bSJed Brown 
217c4762a1bSJed Brown    test:
218c4762a1bSJed Brown       suffix: 32
219c4762a1bSJed Brown       args: -mat_type mpibaij -test_diagonalscale
220c4762a1bSJed Brown 
221c4762a1bSJed Brown    test:
222c4762a1bSJed Brown       suffix: 33
223c4762a1bSJed Brown       nsize: 3
224c4762a1bSJed Brown       args: -mat_type mpiaij -test_diagonalscale
225c4762a1bSJed Brown       filter: grep -v type
226c4762a1bSJed Brown 
227c4762a1bSJed Brown    test:
228c4762a1bSJed Brown       suffix: 34
229c4762a1bSJed Brown       nsize: 3
230c4762a1bSJed Brown       args: -mat_type mpibaij -test_diagonalscale
231c4762a1bSJed Brown 
232c4762a1bSJed Brown    test:
233c4762a1bSJed Brown       suffix: 3_aijcusparse_1
234c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda -test_diagonalscale
235c4762a1bSJed Brown       filter: grep -v type
236c4762a1bSJed Brown       output_file: output/ex5_31.out
237c4762a1bSJed Brown       requires: cuda
238c4762a1bSJed Brown 
239c4762a1bSJed Brown    test:
240c4762a1bSJed Brown       suffix: 3_aijcusparse_2
241c4762a1bSJed Brown       nsize: 3
242c4762a1bSJed Brown       args: -mat_type mpiaijcusparse -vec_type cuda -test_diagonalscale
243c4762a1bSJed Brown       filter: grep -v type
244c4762a1bSJed Brown       output_file: output/ex5_33.out
245c4762a1bSJed Brown       requires: cuda
246c4762a1bSJed Brown 
247c4762a1bSJed Brown    test:
24835990778SJunchao Zhang       suffix: 3_kokkos
24935990778SJunchao Zhang       nsize: 3
25035990778SJunchao Zhang       args: -mat_type mpiaijkokkos -vec_type kokkos -test_diagonalscale
25135990778SJunchao Zhang       filter: grep -v type
25235990778SJunchao Zhang       output_file: output/ex5_33.out
253dcfd994dSJunchao Zhang       requires: kokkos_kernels
25435990778SJunchao Zhang 
25535990778SJunchao Zhang    test:
256c4762a1bSJed Brown       suffix: aijcusparse_1
257c4762a1bSJed Brown       args: -mat_type seqaijcusparse -vec_type cuda -rectA
258c4762a1bSJed Brown       filter: grep -v type
259c4762a1bSJed Brown       output_file: output/ex5_11_A.out
260c4762a1bSJed Brown       requires: cuda
261c4762a1bSJed Brown 
262c4762a1bSJed Brown    test:
263c4762a1bSJed Brown       suffix: aijcusparse_2
264c4762a1bSJed Brown       args: -mat_type seqaijcusparse -vec_type cuda -rectB
265c4762a1bSJed Brown       filter: grep -v type
266c4762a1bSJed Brown       output_file: output/ex5_11_B.out
267c4762a1bSJed Brown       requires: cuda
268c4762a1bSJed Brown 
269c4762a1bSJed Brown    test:
270c4762a1bSJed Brown       suffix: sell_1
271c4762a1bSJed Brown       args: -mat_type sell
272c4762a1bSJed Brown       output_file: output/ex5_41.out
273c4762a1bSJed Brown 
274c4762a1bSJed Brown    test:
275c4762a1bSJed Brown       suffix: sell_2
276c4762a1bSJed Brown       nsize: 3
277c4762a1bSJed Brown       args: -mat_type sell
278c4762a1bSJed Brown       output_file: output/ex5_43.out
279c4762a1bSJed Brown 
280c4762a1bSJed Brown    test:
281c4762a1bSJed Brown       suffix: sell_3
282c4762a1bSJed Brown       args: -mat_type sell -test_diagonalscale
283c4762a1bSJed Brown       output_file: output/ex5_51.out
284c4762a1bSJed Brown 
285c4762a1bSJed Brown    test:
286c4762a1bSJed Brown       suffix: sell_4
287c4762a1bSJed Brown       nsize: 3
288c4762a1bSJed Brown       args: -mat_type sell -test_diagonalscale
289c4762a1bSJed Brown       output_file: output/ex5_53.out
290c4762a1bSJed Brown 
2912d1451d4SHong Zhang    test:
2922d1451d4SHong Zhang       suffix: sell_5
29390d2215bSHong Zhang       nsize: 3
29490d2215bSHong Zhang       args: -mat_type sellcuda -vec_type cuda -test_diagonalscale
29590d2215bSHong Zhang       output_file: output/ex5_55.out
296*8711c661SHong Zhang       requires: cuda !complex
2972d1451d4SHong Zhang 
298c4762a1bSJed Brown TEST*/
299