xref: /petsc/src/mat/tests/ex93.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static char help[] = "Test MatMatMult() and MatPtAP() for AIJ matrices.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown extern PetscErrorCode testPTAPRectangular(void);
6c4762a1bSJed Brown 
79371c9d4SSatish Balay int main(int argc, char **argv) {
8c4762a1bSJed Brown   Mat         A, B, C, D;
9c4762a1bSJed Brown   PetscScalar a[]  = {1., 1., 0., 0., 1., 1., 0., 0., 1.};
10c4762a1bSJed Brown   PetscInt    ij[] = {0, 1, 2};
11c20d7725SJed Brown   PetscReal   fill = 4.0;
12c4762a1bSJed Brown   PetscMPIInt size, rank;
13c20d7725SJed Brown   PetscBool   isequal;
14c20d7725SJed Brown #if defined(PETSC_HAVE_HYPRE)
15c4762a1bSJed Brown   PetscBool test_hypre = PETSC_FALSE;
16c20d7725SJed Brown #endif
17c4762a1bSJed Brown 
18327415f7SBarry Smith   PetscFunctionBeginUser;
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
20c4762a1bSJed Brown #if defined(PETSC_HAVE_HYPRE)
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_hypre", &test_hypre, NULL));
22c4762a1bSJed Brown #endif
239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
25c4762a1bSJed Brown 
269566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
279566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 3, 3));
289566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATAIJ));
299566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
309566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
319566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
32c4762a1bSJed Brown 
33*48a46eb9SPierre Jolivet   if (rank == 0) PetscCall(MatSetValues(A, 3, ij, 3, ij, a, ADD_VALUES));
349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* Test MatMatMult() */
389566063dSJacob Faibussowitsch   PetscCall(MatTranspose(A, MAT_INITIAL_MATRIX, &B));        /* B = A^T */
399566063dSJacob Faibussowitsch   PetscCall(MatMatMult(B, A, MAT_INITIAL_MATRIX, fill, &C)); /* C = B*A */
409566063dSJacob Faibussowitsch   PetscCall(MatMatMult(B, A, MAT_REUSE_MATRIX, fill, &C));   /* recompute C=B*A */
419566063dSJacob Faibussowitsch   PetscCall(MatSetOptionsPrefix(C, "C_"));
429566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(B, A, C, 10, &isequal));
4328b400f6SJacob Faibussowitsch   PetscCheck(isequal, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "MatMatMult: C != B*A");
44c4762a1bSJed Brown 
459566063dSJacob Faibussowitsch   PetscCall(MatMatMult(C, A, MAT_INITIAL_MATRIX, fill, &D)); /* D = C*A = (A^T*A)*A */
469566063dSJacob Faibussowitsch   PetscCall(MatMatMult(C, A, MAT_REUSE_MATRIX, fill, &D));
479566063dSJacob Faibussowitsch   PetscCall(MatMatMultEqual(C, A, D, 10, &isequal));
4828b400f6SJacob Faibussowitsch   PetscCheck(isequal, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "MatMatMult: D != C*A");
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
529566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
53c4762a1bSJed Brown 
54c20d7725SJed Brown   /* Test MatPtAP */
559566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &B));        /* B = A */
569566063dSJacob Faibussowitsch   PetscCall(MatPtAP(A, B, MAT_INITIAL_MATRIX, fill, &C)); /* C = B^T*A*B */
579566063dSJacob Faibussowitsch   PetscCall(MatPtAPMultEqual(A, B, C, 10, &isequal));
5828b400f6SJacob Faibussowitsch   PetscCheck(isequal, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "MatPtAP: C != B^T*A*B");
59c4762a1bSJed Brown 
60c20d7725SJed Brown   /* Repeat MatPtAP to test symbolic/numeric separation for reuse of the symbolic product */
619566063dSJacob Faibussowitsch   PetscCall(MatPtAP(A, B, MAT_REUSE_MATRIX, fill, &C));
629566063dSJacob Faibussowitsch   PetscCall(MatPtAPMultEqual(A, B, C, 10, &isequal));
6328b400f6SJacob Faibussowitsch   PetscCheck(isequal, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "MatPtAP(reuse): C != B^T*A*B");
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
66c1995e1cSHong Zhang 
67c1995e1cSHong Zhang   /* Test MatPtAP with A as a dense matrix */
68c1995e1cSHong Zhang   {
69c1995e1cSHong Zhang     Mat Adense;
709566063dSJacob Faibussowitsch     PetscCall(MatConvert(A, MATDENSE, MAT_INITIAL_MATRIX, &Adense));
719566063dSJacob Faibussowitsch     PetscCall(MatPtAP(Adense, B, MAT_INITIAL_MATRIX, fill, &C));
72c1995e1cSHong Zhang 
739566063dSJacob Faibussowitsch     PetscCall(MatPtAPMultEqual(Adense, B, C, 10, &isequal));
74c1995e1cSHong Zhang     PetscCheck(isequal, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "MatPtAP(reuse): C != B^T*Adense*B");
759566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Adense));
76c1995e1cSHong Zhang   }
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   if (size == 1) {
79c4762a1bSJed Brown     /* A test contributed by Tobias Neckel <neckel@in.tum.de> */
809566063dSJacob Faibussowitsch     PetscCall(testPTAPRectangular());
81c4762a1bSJed Brown 
82c4762a1bSJed Brown     /* test MatMatTransposeMult(): A*B^T */
839566063dSJacob Faibussowitsch     PetscCall(MatMatTransposeMult(A, A, MAT_INITIAL_MATRIX, fill, &D)); /* D = A*A^T */
849566063dSJacob Faibussowitsch     PetscCall(MatScale(A, 2.0));
859566063dSJacob Faibussowitsch     PetscCall(MatMatTransposeMult(A, A, MAT_REUSE_MATRIX, fill, &D));
869566063dSJacob Faibussowitsch     PetscCall(MatMatTransposeMultEqual(A, A, D, 10, &isequal));
8728b400f6SJacob Faibussowitsch     PetscCheck(isequal, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "MatMatTranspose: D != A*A^T");
88c4762a1bSJed Brown   }
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
919566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
929566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&D));
949566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
95b122ec5aSJacob Faibussowitsch   return 0;
96c4762a1bSJed Brown }
97c4762a1bSJed Brown 
98c4762a1bSJed Brown /* a test contributed by Tobias Neckel <neckel@in.tum.de>, 02 Jul 2008 */
999371c9d4SSatish Balay PetscErrorCode testPTAPRectangular(void) {
100c20d7725SJed Brown   const int rows = 3, cols = 5;
101c4762a1bSJed Brown   int       i;
102c20d7725SJed Brown   Mat       A, P, C;
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   PetscFunctionBegin;
105c4762a1bSJed Brown   /* set up A  */
1069566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, rows, rows, 1, NULL, &A));
107*48a46eb9SPierre Jolivet   for (i = 0; i < rows; i++) PetscCall(MatSetValue(A, i, i, 1.0, INSERT_VALUES));
1089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1099566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   /* set up P */
1129566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, rows, cols, 5, NULL, &P));
1139566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P, 0, 0, 1.0, INSERT_VALUES));
1149566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P, 0, 1, 2.0, INSERT_VALUES));
1159566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P, 0, 2, 0.0, INSERT_VALUES));
116c4762a1bSJed Brown 
1179566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P, 0, 3, -1.0, INSERT_VALUES));
118c4762a1bSJed Brown 
1199566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P, 1, 0, 0.0, INSERT_VALUES));
1209566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P, 1, 1, -1.0, INSERT_VALUES));
1219566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P, 1, 2, 1.0, INSERT_VALUES));
122c4762a1bSJed Brown 
1239566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P, 2, 0, 3.0, INSERT_VALUES));
1249566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P, 2, 1, 0.0, INSERT_VALUES));
1259566063dSJacob Faibussowitsch   PetscCall(MatSetValue(P, 2, 2, -3.0, INSERT_VALUES));
126c4762a1bSJed Brown 
1279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(P, MAT_FINAL_ASSEMBLY));
1289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(P, MAT_FINAL_ASSEMBLY));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* compute C */
1319566063dSJacob Faibussowitsch   PetscCall(MatPtAP(A, P, MAT_INITIAL_MATRIX, 1.0, &C));
132c4762a1bSJed Brown 
1339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
1349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* compare results */
137c4762a1bSJed Brown   /*
138c4762a1bSJed Brown   printf("C:\n");
1399566063dSJacob Faibussowitsch   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   blitz::Array<double,2> actualC(cols, cols);
142c4762a1bSJed Brown   actualC = 0.0;
143c4762a1bSJed Brown   for (int i=0; i<cols; i++) {
144c4762a1bSJed Brown     for (int j=0; j<cols; j++) {
1459566063dSJacob Faibussowitsch       PetscCall(MatGetValues(C, 1, &i, 1, &j, &actualC(i,j)));
146c4762a1bSJed Brown     }
147c4762a1bSJed Brown   }
148c4762a1bSJed Brown   blitz::Array<double,2> expectedC(cols, cols);
149c4762a1bSJed Brown   expectedC = 0.0;
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   expectedC(0,0) = 10.0;
152c4762a1bSJed Brown   expectedC(0,1) =  2.0;
153c4762a1bSJed Brown   expectedC(0,2) = -9.0;
154c4762a1bSJed Brown   expectedC(0,3) = -1.0;
155c4762a1bSJed Brown   expectedC(1,0) =  2.0;
156c4762a1bSJed Brown   expectedC(1,1) =  5.0;
157c4762a1bSJed Brown   expectedC(1,2) = -1.0;
158c4762a1bSJed Brown   expectedC(1,3) = -2.0;
159c4762a1bSJed Brown   expectedC(2,0) = -9.0;
160c4762a1bSJed Brown   expectedC(2,1) = -1.0;
161c4762a1bSJed Brown   expectedC(2,2) = 10.0;
162c4762a1bSJed Brown   expectedC(2,3) =  0.0;
163c4762a1bSJed Brown   expectedC(3,0) = -1.0;
164c4762a1bSJed Brown   expectedC(3,1) = -2.0;
165c4762a1bSJed Brown   expectedC(3,2) =  0.0;
166c4762a1bSJed Brown   expectedC(3,3) =  1.0;
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   int check = areBlitzArrays2NumericallyEqual(actualC,expectedC);
169c4762a1bSJed Brown   validateEqualsWithParams3(check, -1 , "testPTAPRectangular()", check, actualC(check), expectedC(check));
170c4762a1bSJed Brown   */
171c4762a1bSJed Brown 
1729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&P));
1749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
175c4762a1bSJed Brown   PetscFunctionReturn(0);
176c4762a1bSJed Brown }
177c4762a1bSJed Brown 
178c4762a1bSJed Brown /*TEST
179c4762a1bSJed Brown 
180c4762a1bSJed Brown    test:
181c4762a1bSJed Brown 
182c4762a1bSJed Brown    test:
183c4762a1bSJed Brown       suffix: 2
184c4762a1bSJed Brown       nsize: 2
185c20d7725SJed Brown       args: -matmatmult_via nonscalable
186c20d7725SJed Brown       output_file: output/ex93_1.out
187c4762a1bSJed Brown 
188c4762a1bSJed Brown    test:
189c4762a1bSJed Brown       suffix: 3
190c4762a1bSJed Brown       nsize: 2
191c20d7725SJed Brown       output_file: output/ex93_1.out
192c4762a1bSJed Brown 
193c4762a1bSJed Brown    test:
194c4762a1bSJed Brown       suffix: 4
195c4762a1bSJed Brown       nsize: 2
196c20d7725SJed Brown       args: -matptap_via scalable
197c20d7725SJed Brown       output_file: output/ex93_1.out
198c4762a1bSJed Brown 
199c4762a1bSJed Brown    test:
200c4762a1bSJed Brown       suffix: btheap
201c20d7725SJed Brown       args: -matmatmult_via btheap -matmattransmult_via color
202c4762a1bSJed Brown       output_file: output/ex93_1.out
203c4762a1bSJed Brown 
204c4762a1bSJed Brown    test:
205c4762a1bSJed Brown       suffix: heap
206c20d7725SJed Brown       args: -matmatmult_via heap
207c4762a1bSJed Brown       output_file: output/ex93_1.out
208c4762a1bSJed Brown 
209c4762a1bSJed Brown    #HYPRE PtAP is broken for complex numbers
210c4762a1bSJed Brown    test:
211c4762a1bSJed Brown       suffix: hypre
212c4762a1bSJed Brown       nsize: 3
213c4762a1bSJed Brown       requires: hypre !complex
214c20d7725SJed Brown       args: -matmatmult_via hypre -matptap_via hypre -test_hypre
215c20d7725SJed Brown       output_file: output/ex93_hypre.out
216c4762a1bSJed Brown 
217c4762a1bSJed Brown    test:
218c4762a1bSJed Brown       suffix: llcondensed
219c20d7725SJed Brown       args: -matmatmult_via llcondensed
220c4762a1bSJed Brown       output_file: output/ex93_1.out
221c4762a1bSJed Brown 
222c4762a1bSJed Brown    test:
223c4762a1bSJed Brown       suffix: scalable
224c20d7725SJed Brown       args: -matmatmult_via scalable
225c4762a1bSJed Brown       output_file: output/ex93_1.out
226c4762a1bSJed Brown 
227c4762a1bSJed Brown    test:
228c4762a1bSJed Brown       suffix: scalable_fast
229c20d7725SJed Brown       args: -matmatmult_via scalable_fast
230c4762a1bSJed Brown       output_file: output/ex93_1.out
231c4762a1bSJed Brown 
232c4762a1bSJed Brown TEST*/
233