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