xref: /petsc/src/mat/tests/ex93.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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 
7c4762a1bSJed Brown int main(int argc,char **argv)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   Mat            A,B,C,D;
10c4762a1bSJed Brown   PetscScalar    a[] ={1.,1.,0.,0.,1.,1.,0.,0.,1.};
11c4762a1bSJed Brown   PetscInt       ij[]={0,1,2};
12c4762a1bSJed Brown   PetscErrorCode ierr;
13c20d7725SJed Brown   PetscReal      fill=4.0;
14c4762a1bSJed Brown   PetscMPIInt    size,rank;
15c20d7725SJed Brown   PetscBool      isequal;
16c20d7725SJed Brown #if defined(PETSC_HAVE_HYPRE)
17c4762a1bSJed Brown   PetscBool      test_hypre=PETSC_FALSE;
18c20d7725SJed Brown #endif
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
21c4762a1bSJed Brown #if defined(PETSC_HAVE_HYPRE)
225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_hypre",&test_hypre,NULL));
23c4762a1bSJed Brown #endif
245f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
255f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
26c4762a1bSJed Brown 
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,3,3));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATAIJ));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE));
33c4762a1bSJed Brown 
34dd400576SPatrick Sanan   if (rank == 0) {
355f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,3,ij,3,ij,a,ADD_VALUES));
36c4762a1bSJed Brown   }
375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   /* Test MatMatMult() */
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&B));      /* B = A^T */
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C)); /* C = B*A */
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C));   /* recompute C=B*A */
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOptionsPrefix(C,"C_"));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(B,A,C,10,&isequal));
46*28b400f6SJacob Faibussowitsch   PetscCheck(isequal,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatMatMult: C != B*A");
47c4762a1bSJed Brown 
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D)); /* D = C*A = (A^T*A)*A */
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMult(C,A,MAT_REUSE_MATRIX,fill,&D));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatMultEqual(C,A,D,10,&isequal));
51*28b400f6SJacob Faibussowitsch   PetscCheck(isequal,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatMatMult: D != C*A");
52c4762a1bSJed Brown 
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
56c4762a1bSJed Brown 
57c20d7725SJed Brown   /* Test MatPtAP */
585f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&B));      /* B = A */
595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPtAP(A,B,MAT_INITIAL_MATRIX,fill,&C)); /* C = B^T*A*B */
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPtAPMultEqual(A,B,C,10,&isequal));
61*28b400f6SJacob Faibussowitsch   PetscCheck(isequal,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatPtAP: C != B^T*A*B");
62c4762a1bSJed Brown 
63c20d7725SJed Brown   /* Repeat MatPtAP to test symbolic/numeric separation for reuse of the symbolic product */
645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPtAP(A,B,MAT_REUSE_MATRIX,fill,&C));
655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPtAPMultEqual(A,B,C,10,&isequal));
66*28b400f6SJacob Faibussowitsch   PetscCheck(isequal,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatPtAP(reuse): C != B^T*A*B");
67c4762a1bSJed Brown 
685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
69c1995e1cSHong Zhang 
70c1995e1cSHong Zhang   /* Test MatPtAP with A as a dense matrix */
71c1995e1cSHong Zhang   {
72c1995e1cSHong Zhang     Mat Adense;
735f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense));
745f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPtAP(Adense,B,MAT_INITIAL_MATRIX,fill,&C));
75c1995e1cSHong Zhang 
765f80ce2aSJacob Faibussowitsch     CHKERRQ(MatPtAPMultEqual(Adense,B,C,10,&isequal));
77c1995e1cSHong Zhang     PetscCheck(isequal,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatPtAP(reuse): C != B^T*Adense*B");
785f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&Adense));
79c1995e1cSHong Zhang   }
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   if (size == 1) {
82c4762a1bSJed Brown     /* A test contributed by Tobias Neckel <neckel@in.tum.de> */
835f80ce2aSJacob Faibussowitsch     CHKERRQ(testPTAPRectangular());
84c4762a1bSJed Brown 
85c4762a1bSJed Brown     /* test MatMatTransposeMult(): A*B^T */
865f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatTransposeMult(A,A,MAT_INITIAL_MATRIX,fill,&D)); /* D = A*A^T */
875f80ce2aSJacob Faibussowitsch     CHKERRQ(MatScale(A,2.0));
885f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatTransposeMult(A,A,MAT_REUSE_MATRIX,fill,&D));
895f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatTransposeMultEqual(A,A,D,10,&isequal));
90*28b400f6SJacob Faibussowitsch     PetscCheck(isequal,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatMatTranspose: D != A*A^T");
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown 
935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&D));
97c4762a1bSJed Brown   ierr = PetscFinalize();
98c4762a1bSJed Brown   return ierr;
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown /* a test contributed by Tobias Neckel <neckel@in.tum.de>, 02 Jul 2008 */
102c4762a1bSJed Brown PetscErrorCode testPTAPRectangular(void)
103c4762a1bSJed Brown {
104c20d7725SJed Brown   const int      rows = 3,cols = 5;
105c4762a1bSJed Brown   int            i;
106c20d7725SJed Brown   Mat            A,P,C;
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   PetscFunctionBegin;
109c4762a1bSJed Brown   /* set up A  */
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD, rows, rows,1, NULL, &A));
111c4762a1bSJed Brown   for (i=0; i<rows; i++) {
1125f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValue(A, i, i, 1.0, INSERT_VALUES));
113c4762a1bSJed Brown   }
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /* set up P */
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD, rows, cols,5, NULL, &P));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(P, 0, 0,  1.0, INSERT_VALUES));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(P, 0, 1,  2.0, INSERT_VALUES));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(P, 0, 2,  0.0, INSERT_VALUES));
122c4762a1bSJed Brown 
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(P, 0, 3, -1.0, INSERT_VALUES));
124c4762a1bSJed Brown 
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(P, 1, 0,  0.0, INSERT_VALUES));
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(P, 1, 1, -1.0, INSERT_VALUES));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(P, 1, 2,  1.0, INSERT_VALUES));
128c4762a1bSJed Brown 
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(P, 2, 0,  3.0, INSERT_VALUES));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(P, 2, 1,  0.0, INSERT_VALUES));
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(P, 2, 2, -3.0, INSERT_VALUES));
132c4762a1bSJed Brown 
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY));
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY));
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* compute C */
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatPtAP(A, P, MAT_INITIAL_MATRIX, 1.0, &C));
138c4762a1bSJed Brown 
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
1405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* compare results */
143c4762a1bSJed Brown   /*
144c4762a1bSJed Brown   printf("C:\n");
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   blitz::Array<double,2> actualC(cols, cols);
148c4762a1bSJed Brown   actualC = 0.0;
149c4762a1bSJed Brown   for (int i=0; i<cols; i++) {
150c4762a1bSJed Brown     for (int j=0; j<cols; j++) {
1515f80ce2aSJacob Faibussowitsch       CHKERRQ(MatGetValues(C, 1, &i, 1, &j, &actualC(i,j)));
152c4762a1bSJed Brown     }
153c4762a1bSJed Brown   }
154c4762a1bSJed Brown   blitz::Array<double,2> expectedC(cols, cols);
155c4762a1bSJed Brown   expectedC = 0.0;
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   expectedC(0,0) = 10.0;
158c4762a1bSJed Brown   expectedC(0,1) =  2.0;
159c4762a1bSJed Brown   expectedC(0,2) = -9.0;
160c4762a1bSJed Brown   expectedC(0,3) = -1.0;
161c4762a1bSJed Brown   expectedC(1,0) =  2.0;
162c4762a1bSJed Brown   expectedC(1,1) =  5.0;
163c4762a1bSJed Brown   expectedC(1,2) = -1.0;
164c4762a1bSJed Brown   expectedC(1,3) = -2.0;
165c4762a1bSJed Brown   expectedC(2,0) = -9.0;
166c4762a1bSJed Brown   expectedC(2,1) = -1.0;
167c4762a1bSJed Brown   expectedC(2,2) = 10.0;
168c4762a1bSJed Brown   expectedC(2,3) =  0.0;
169c4762a1bSJed Brown   expectedC(3,0) = -1.0;
170c4762a1bSJed Brown   expectedC(3,1) = -2.0;
171c4762a1bSJed Brown   expectedC(3,2) =  0.0;
172c4762a1bSJed Brown   expectedC(3,3) =  1.0;
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   int check = areBlitzArrays2NumericallyEqual(actualC,expectedC);
175c4762a1bSJed Brown   validateEqualsWithParams3(check, -1 , "testPTAPRectangular()", check, actualC(check), expectedC(check));
176c4762a1bSJed Brown   */
177c4762a1bSJed Brown 
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&P));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
181c4762a1bSJed Brown   PetscFunctionReturn(0);
182c4762a1bSJed Brown }
183c4762a1bSJed Brown 
184c4762a1bSJed Brown /*TEST
185c4762a1bSJed Brown 
186c4762a1bSJed Brown    test:
187c4762a1bSJed Brown 
188c4762a1bSJed Brown    test:
189c4762a1bSJed Brown       suffix: 2
190c4762a1bSJed Brown       nsize: 2
191c20d7725SJed Brown       args: -matmatmult_via nonscalable
192c20d7725SJed Brown       output_file: output/ex93_1.out
193c4762a1bSJed Brown 
194c4762a1bSJed Brown    test:
195c4762a1bSJed Brown       suffix: 3
196c4762a1bSJed Brown       nsize: 2
197c20d7725SJed Brown       output_file: output/ex93_1.out
198c4762a1bSJed Brown 
199c4762a1bSJed Brown    test:
200c4762a1bSJed Brown       suffix: 4
201c4762a1bSJed Brown       nsize: 2
202c20d7725SJed Brown       args: -matptap_via scalable
203c20d7725SJed Brown       output_file: output/ex93_1.out
204c4762a1bSJed Brown 
205c4762a1bSJed Brown    test:
206c4762a1bSJed Brown       suffix: btheap
207c20d7725SJed Brown       args: -matmatmult_via btheap -matmattransmult_via color
208c4762a1bSJed Brown       output_file: output/ex93_1.out
209c4762a1bSJed Brown 
210c4762a1bSJed Brown    test:
211c4762a1bSJed Brown       suffix: heap
212c20d7725SJed Brown       args: -matmatmult_via heap
213c4762a1bSJed Brown       output_file: output/ex93_1.out
214c4762a1bSJed Brown 
215c4762a1bSJed Brown    #HYPRE PtAP is broken for complex numbers
216c4762a1bSJed Brown    test:
217c4762a1bSJed Brown       suffix: hypre
218c4762a1bSJed Brown       nsize: 3
219c4762a1bSJed Brown       requires: hypre !complex
220c20d7725SJed Brown       args: -matmatmult_via hypre -matptap_via hypre -test_hypre
221c20d7725SJed Brown       output_file: output/ex93_hypre.out
222c4762a1bSJed Brown 
223c4762a1bSJed Brown    test:
224c4762a1bSJed Brown       suffix: llcondensed
225c20d7725SJed Brown       args: -matmatmult_via llcondensed
226c4762a1bSJed Brown       output_file: output/ex93_1.out
227c4762a1bSJed Brown 
228c4762a1bSJed Brown    test:
229c4762a1bSJed Brown       suffix: scalable
230c20d7725SJed Brown       args: -matmatmult_via scalable
231c4762a1bSJed Brown       output_file: output/ex93_1.out
232c4762a1bSJed Brown 
233c4762a1bSJed Brown    test:
234c4762a1bSJed Brown       suffix: scalable_fast
235c20d7725SJed Brown       args: -matmatmult_via scalable_fast
236c4762a1bSJed Brown       output_file: output/ex93_1.out
237c4762a1bSJed Brown 
238c4762a1bSJed Brown TEST*/
239