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}; 12c20d7725SJed Brown PetscReal fill=4.0; 13c4762a1bSJed Brown PetscMPIInt size,rank; 14c20d7725SJed Brown PetscBool isequal; 15c20d7725SJed Brown #if defined(PETSC_HAVE_HYPRE) 16c4762a1bSJed Brown PetscBool test_hypre=PETSC_FALSE; 17c20d7725SJed Brown #endif 18c4762a1bSJed Brown 19*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 20c4762a1bSJed Brown #if defined(PETSC_HAVE_HYPRE) 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_hypre",&test_hypre,NULL)); 22c4762a1bSJed Brown #endif 235f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 245f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 25c4762a1bSJed Brown 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,3,3)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATAIJ)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE)); 32c4762a1bSJed Brown 33dd400576SPatrick Sanan if (rank == 0) { 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,3,ij,3,ij,a,ADD_VALUES)); 35c4762a1bSJed Brown } 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Test MatMatMult() */ 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&B)); /* B = A^T */ 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C)); /* C = B*A */ 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C)); /* recompute C=B*A */ 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(C,"C_")); 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(B,A,C,10,&isequal)); 4528b400f6SJacob Faibussowitsch PetscCheck(isequal,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatMatMult: C != B*A"); 46c4762a1bSJed Brown 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D)); /* D = C*A = (A^T*A)*A */ 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(C,A,MAT_REUSE_MATRIX,fill,&D)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMultEqual(C,A,D,10,&isequal)); 5028b400f6SJacob Faibussowitsch PetscCheck(isequal,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatMatMult: D != C*A"); 51c4762a1bSJed Brown 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 55c4762a1bSJed Brown 56c20d7725SJed Brown /* Test MatPtAP */ 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&B)); /* B = A */ 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,B,MAT_INITIAL_MATRIX,fill,&C)); /* C = B^T*A*B */ 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAPMultEqual(A,B,C,10,&isequal)); 6028b400f6SJacob Faibussowitsch PetscCheck(isequal,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatPtAP: C != B^T*A*B"); 61c4762a1bSJed Brown 62c20d7725SJed Brown /* Repeat MatPtAP to test symbolic/numeric separation for reuse of the symbolic product */ 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A,B,MAT_REUSE_MATRIX,fill,&C)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAPMultEqual(A,B,C,10,&isequal)); 6528b400f6SJacob Faibussowitsch PetscCheck(isequal,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatPtAP(reuse): C != B^T*A*B"); 66c4762a1bSJed Brown 675f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 68c1995e1cSHong Zhang 69c1995e1cSHong Zhang /* Test MatPtAP with A as a dense matrix */ 70c1995e1cSHong Zhang { 71c1995e1cSHong Zhang Mat Adense; 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATDENSE,MAT_INITIAL_MATRIX,&Adense)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(Adense,B,MAT_INITIAL_MATRIX,fill,&C)); 74c1995e1cSHong Zhang 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAPMultEqual(Adense,B,C,10,&isequal)); 76c1995e1cSHong Zhang PetscCheck(isequal,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatPtAP(reuse): C != B^T*Adense*B"); 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Adense)); 78c1995e1cSHong Zhang } 79c4762a1bSJed Brown 80c4762a1bSJed Brown if (size == 1) { 81c4762a1bSJed Brown /* A test contributed by Tobias Neckel <neckel@in.tum.de> */ 825f80ce2aSJacob Faibussowitsch CHKERRQ(testPTAPRectangular()); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* test MatMatTransposeMult(): A*B^T */ 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(A,A,MAT_INITIAL_MATRIX,fill,&D)); /* D = A*A^T */ 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(A,2.0)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(A,A,MAT_REUSE_MATRIX,fill,&D)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMultEqual(A,A,D,10,&isequal)); 8928b400f6SJacob Faibussowitsch PetscCheck(isequal,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"MatMatTranspose: D != A*A^T"); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&D)); 96*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 97*b122ec5aSJacob Faibussowitsch return 0; 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* a test contributed by Tobias Neckel <neckel@in.tum.de>, 02 Jul 2008 */ 101c4762a1bSJed Brown PetscErrorCode testPTAPRectangular(void) 102c4762a1bSJed Brown { 103c20d7725SJed Brown const int rows = 3,cols = 5; 104c4762a1bSJed Brown int i; 105c20d7725SJed Brown Mat A,P,C; 106c4762a1bSJed Brown 107c4762a1bSJed Brown PetscFunctionBegin; 108c4762a1bSJed Brown /* set up A */ 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD, rows, rows,1, NULL, &A)); 110c4762a1bSJed Brown for (i=0; i<rows; i++) { 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(A, i, i, 1.0, INSERT_VALUES)); 112c4762a1bSJed Brown } 1135f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* set up P */ 1175f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD, rows, cols,5, NULL, &P)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(P, 0, 0, 1.0, INSERT_VALUES)); 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(P, 0, 1, 2.0, INSERT_VALUES)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(P, 0, 2, 0.0, INSERT_VALUES)); 121c4762a1bSJed Brown 1225f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(P, 0, 3, -1.0, INSERT_VALUES)); 123c4762a1bSJed Brown 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(P, 1, 0, 0.0, INSERT_VALUES)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(P, 1, 1, -1.0, INSERT_VALUES)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(P, 1, 2, 1.0, INSERT_VALUES)); 127c4762a1bSJed Brown 1285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(P, 2, 0, 3.0, INSERT_VALUES)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(P, 2, 1, 0.0, INSERT_VALUES)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValue(P, 2, 2, -3.0, INSERT_VALUES)); 131c4762a1bSJed Brown 1325f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY)); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* compute C */ 1365f80ce2aSJacob Faibussowitsch CHKERRQ(MatPtAP(A, P, MAT_INITIAL_MATRIX, 1.0, &C)); 137c4762a1bSJed Brown 1385f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* compare results */ 142c4762a1bSJed Brown /* 143c4762a1bSJed Brown printf("C:\n"); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 145c4762a1bSJed Brown 146c4762a1bSJed Brown blitz::Array<double,2> actualC(cols, cols); 147c4762a1bSJed Brown actualC = 0.0; 148c4762a1bSJed Brown for (int i=0; i<cols; i++) { 149c4762a1bSJed Brown for (int j=0; j<cols; j++) { 1505f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetValues(C, 1, &i, 1, &j, &actualC(i,j))); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown } 153c4762a1bSJed Brown blitz::Array<double,2> expectedC(cols, cols); 154c4762a1bSJed Brown expectedC = 0.0; 155c4762a1bSJed Brown 156c4762a1bSJed Brown expectedC(0,0) = 10.0; 157c4762a1bSJed Brown expectedC(0,1) = 2.0; 158c4762a1bSJed Brown expectedC(0,2) = -9.0; 159c4762a1bSJed Brown expectedC(0,3) = -1.0; 160c4762a1bSJed Brown expectedC(1,0) = 2.0; 161c4762a1bSJed Brown expectedC(1,1) = 5.0; 162c4762a1bSJed Brown expectedC(1,2) = -1.0; 163c4762a1bSJed Brown expectedC(1,3) = -2.0; 164c4762a1bSJed Brown expectedC(2,0) = -9.0; 165c4762a1bSJed Brown expectedC(2,1) = -1.0; 166c4762a1bSJed Brown expectedC(2,2) = 10.0; 167c4762a1bSJed Brown expectedC(2,3) = 0.0; 168c4762a1bSJed Brown expectedC(3,0) = -1.0; 169c4762a1bSJed Brown expectedC(3,1) = -2.0; 170c4762a1bSJed Brown expectedC(3,2) = 0.0; 171c4762a1bSJed Brown expectedC(3,3) = 1.0; 172c4762a1bSJed Brown 173c4762a1bSJed Brown int check = areBlitzArrays2NumericallyEqual(actualC,expectedC); 174c4762a1bSJed Brown validateEqualsWithParams3(check, -1 , "testPTAPRectangular()", check, actualC(check), expectedC(check)); 175c4762a1bSJed Brown */ 176c4762a1bSJed Brown 1775f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&P)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 180c4762a1bSJed Brown PetscFunctionReturn(0); 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183c4762a1bSJed Brown /*TEST 184c4762a1bSJed Brown 185c4762a1bSJed Brown test: 186c4762a1bSJed Brown 187c4762a1bSJed Brown test: 188c4762a1bSJed Brown suffix: 2 189c4762a1bSJed Brown nsize: 2 190c20d7725SJed Brown args: -matmatmult_via nonscalable 191c20d7725SJed Brown output_file: output/ex93_1.out 192c4762a1bSJed Brown 193c4762a1bSJed Brown test: 194c4762a1bSJed Brown suffix: 3 195c4762a1bSJed Brown nsize: 2 196c20d7725SJed Brown output_file: output/ex93_1.out 197c4762a1bSJed Brown 198c4762a1bSJed Brown test: 199c4762a1bSJed Brown suffix: 4 200c4762a1bSJed Brown nsize: 2 201c20d7725SJed Brown args: -matptap_via scalable 202c20d7725SJed Brown output_file: output/ex93_1.out 203c4762a1bSJed Brown 204c4762a1bSJed Brown test: 205c4762a1bSJed Brown suffix: btheap 206c20d7725SJed Brown args: -matmatmult_via btheap -matmattransmult_via color 207c4762a1bSJed Brown output_file: output/ex93_1.out 208c4762a1bSJed Brown 209c4762a1bSJed Brown test: 210c4762a1bSJed Brown suffix: heap 211c20d7725SJed Brown args: -matmatmult_via heap 212c4762a1bSJed Brown output_file: output/ex93_1.out 213c4762a1bSJed Brown 214c4762a1bSJed Brown #HYPRE PtAP is broken for complex numbers 215c4762a1bSJed Brown test: 216c4762a1bSJed Brown suffix: hypre 217c4762a1bSJed Brown nsize: 3 218c4762a1bSJed Brown requires: hypre !complex 219c20d7725SJed Brown args: -matmatmult_via hypre -matptap_via hypre -test_hypre 220c20d7725SJed Brown output_file: output/ex93_hypre.out 221c4762a1bSJed Brown 222c4762a1bSJed Brown test: 223c4762a1bSJed Brown suffix: llcondensed 224c20d7725SJed Brown args: -matmatmult_via llcondensed 225c4762a1bSJed Brown output_file: output/ex93_1.out 226c4762a1bSJed Brown 227c4762a1bSJed Brown test: 228c4762a1bSJed Brown suffix: scalable 229c20d7725SJed Brown args: -matmatmult_via scalable 230c4762a1bSJed Brown output_file: output/ex93_1.out 231c4762a1bSJed Brown 232c4762a1bSJed Brown test: 233c4762a1bSJed Brown suffix: scalable_fast 234c20d7725SJed Brown args: -matmatmult_via scalable_fast 235c4762a1bSJed Brown output_file: output/ex93_1.out 236c4762a1bSJed Brown 237c4762a1bSJed Brown TEST*/ 238