1*c4762a1bSJed Brown static char help[] = "Test MatMatMult() and MatPtAP() for AIJ matrices.\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscmat.h> 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown extern PetscErrorCode testPTAPRectangular(void); 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown int main(int argc,char **argv) 8*c4762a1bSJed Brown { 9*c4762a1bSJed Brown Mat A,B,C,D; 10*c4762a1bSJed Brown PetscScalar a[] ={1.,1.,0.,0.,1.,1.,0.,0.,1.}; 11*c4762a1bSJed Brown PetscInt ij[]={0,1,2}; 12*c4762a1bSJed Brown PetscScalar none=-1.; 13*c4762a1bSJed Brown PetscErrorCode ierr; 14*c4762a1bSJed Brown PetscReal fill=4; 15*c4762a1bSJed Brown PetscReal norm; 16*c4762a1bSJed Brown PetscMPIInt size,rank; 17*c4762a1bSJed Brown PetscBool test_hypre=PETSC_FALSE; 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 20*c4762a1bSJed Brown #if defined(PETSC_HAVE_HYPRE) 21*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_hypre",&test_hypre,NULL);CHKERRQ(ierr); 22*c4762a1bSJed Brown #endif 23*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,3,3);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = MatSetOption(A,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown if (!rank) { 34*c4762a1bSJed Brown ierr = MatSetValues(A,3,ij,3,ij,a,ADD_VALUES);CHKERRQ(ierr); 35*c4762a1bSJed Brown } 36*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 38*c4762a1bSJed Brown ierr = MatSetOptionsPrefix(A,"A_");CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n");CHKERRQ(ierr); 41*c4762a1bSJed Brown 42*c4762a1bSJed Brown /* Test MatMatMult() */ 43*c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); /* B = A^T */ 44*c4762a1bSJed Brown ierr = MatSetOptionsPrefix(B,"B_");CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); /* C = B*A */ 46*c4762a1bSJed Brown ierr = MatMatMultNumeric(B,A,C);CHKERRQ(ierr); /* recompute C=B*A */ 47*c4762a1bSJed Brown ierr = MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C);CHKERRQ(ierr); /* recompute C=B*A */ 48*c4762a1bSJed Brown ierr = MatSetOptionsPrefix(C,"C_");CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 50*c4762a1bSJed Brown if (!rank) {ierr = PetscPrintf(PETSC_COMM_SELF,"\n");CHKERRQ(ierr);} 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown ierr = MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = MatMatMultNumeric(C,A,D);CHKERRQ(ierr); /* D = C*A = (A^T*A)*A */ 54*c4762a1bSJed Brown ierr = MatSetOptionsPrefix(D,"D_");CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = MatView(D,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 56*c4762a1bSJed Brown if (!rank) {ierr = PetscPrintf(PETSC_COMM_SELF,"\n");CHKERRQ(ierr);} 57*c4762a1bSJed Brown 58*c4762a1bSJed Brown /* Repeat the numeric product to test reuse of the previous symbolic product */ 59*c4762a1bSJed Brown ierr = MatMatMultNumeric(C,A,D);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = MatView(D,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 61*c4762a1bSJed Brown if (!rank) {ierr = PetscPrintf(PETSC_COMM_SELF,"\n");CHKERRQ(ierr);} 62*c4762a1bSJed Brown 63*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown /* Test PtAP routine. */ 67*c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); /* B = A */ 68*c4762a1bSJed Brown ierr = MatPtAP(A,B,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); /* C = B^T*A*B */ 69*c4762a1bSJed Brown if (!test_hypre) { 70*c4762a1bSJed Brown ierr = MatAXPY(D,none,C,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = MatNorm(D,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 72*c4762a1bSJed Brown if (norm > 1.e-15 && !rank) { 73*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Error in MatPtAP: %g\n",norm);CHKERRQ(ierr); 74*c4762a1bSJed Brown } 75*c4762a1bSJed Brown } 76*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 78*c4762a1bSJed Brown 79*c4762a1bSJed Brown /* Repeat PtAP to test symbolic/numeric separation for reuse of the symbolic product */ 80*c4762a1bSJed Brown ierr = MatPtAP(A,B,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = MatSetOptionsPrefix(C,"C=BtAB_");CHKERRQ(ierr); 82*c4762a1bSJed Brown if (!test_hypre) { 83*c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n");CHKERRQ(ierr); 85*c4762a1bSJed Brown } 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown if (!test_hypre) { 88*c4762a1bSJed Brown ierr = MatPtAPSymbolic(A,B,fill,&D);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = MatPtAPNumeric(A,B,D);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = MatSetOptionsPrefix(D,"D=BtAB_");CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = MatView(D,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n");CHKERRQ(ierr); 93*c4762a1bSJed Brown 94*c4762a1bSJed Brown /* Repeat numeric product to test reuse of the previous symbolic product */ 95*c4762a1bSJed Brown ierr = MatPtAPNumeric(A,B,D);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = MatAXPY(D,none,C,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = MatNorm(D,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 98*c4762a1bSJed Brown if (norm > 1.e-15) { 99*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Error in symbolic/numeric MatPtAP: %g\n",norm);CHKERRQ(ierr); 100*c4762a1bSJed Brown } 101*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 102*c4762a1bSJed Brown } 103*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown if (size == 1) { 107*c4762a1bSJed Brown /* A test contributed by Tobias Neckel <neckel@in.tum.de> */ 108*c4762a1bSJed Brown ierr = testPTAPRectangular();CHKERRQ(ierr); 109*c4762a1bSJed Brown 110*c4762a1bSJed Brown /* test MatMatTransposeMult(): A*B^T */ 111*c4762a1bSJed Brown ierr = MatMatTransposeMult(A,A,MAT_INITIAL_MATRIX,fill,&D);CHKERRQ(ierr); /* D = A*A^T */ 112*c4762a1bSJed Brown ierr = MatSetOptionsPrefix(D,"D=A*A^T_");CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = MatView(D,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n");CHKERRQ(ierr); 115*c4762a1bSJed Brown 116*c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); /* B = A^T */ 117*c4762a1bSJed Brown ierr = MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C);CHKERRQ(ierr); /* C=A*B */ 118*c4762a1bSJed Brown ierr = MatSetOptionsPrefix(C,"D=A*B=A*A^T_");CHKERRQ(ierr); 119*c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n");CHKERRQ(ierr); 121*c4762a1bSJed Brown } 122*c4762a1bSJed Brown 123*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 124*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = MatDestroy(&D);CHKERRQ(ierr); 127*c4762a1bSJed Brown ierr = PetscFinalize(); 128*c4762a1bSJed Brown return ierr; 129*c4762a1bSJed Brown } 130*c4762a1bSJed Brown 131*c4762a1bSJed Brown /* a test contributed by Tobias Neckel <neckel@in.tum.de>, 02 Jul 2008 */ 132*c4762a1bSJed Brown PetscErrorCode testPTAPRectangular(void) 133*c4762a1bSJed Brown { 134*c4762a1bSJed Brown 135*c4762a1bSJed Brown const int rows = 3; 136*c4762a1bSJed Brown const int cols = 5; 137*c4762a1bSJed Brown PetscErrorCode ierr; 138*c4762a1bSJed Brown int i; 139*c4762a1bSJed Brown Mat A; 140*c4762a1bSJed Brown Mat P; 141*c4762a1bSJed Brown Mat C; 142*c4762a1bSJed Brown 143*c4762a1bSJed Brown PetscFunctionBegin; 144*c4762a1bSJed Brown /* set up A */ 145*c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD, rows, rows,1, NULL, &A);CHKERRQ(ierr); 146*c4762a1bSJed Brown for (i=0; i<rows; i++) { 147*c4762a1bSJed Brown ierr = MatSetValue(A, i, i, 1.0, INSERT_VALUES);CHKERRQ(ierr); 148*c4762a1bSJed Brown } 149*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 150*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 151*c4762a1bSJed Brown 152*c4762a1bSJed Brown /* set up P */ 153*c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD, rows, cols,5, NULL, &P);CHKERRQ(ierr); 154*c4762a1bSJed Brown ierr = MatSetValue(P, 0, 0, 1.0, INSERT_VALUES);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = MatSetValue(P, 0, 1, 2.0, INSERT_VALUES);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = MatSetValue(P, 0, 2, 0.0, INSERT_VALUES);CHKERRQ(ierr); 157*c4762a1bSJed Brown 158*c4762a1bSJed Brown ierr = MatSetValue(P, 0, 3, -1.0, INSERT_VALUES);CHKERRQ(ierr); 159*c4762a1bSJed Brown 160*c4762a1bSJed Brown ierr = MatSetValue(P, 1, 0, 0.0, INSERT_VALUES);CHKERRQ(ierr); 161*c4762a1bSJed Brown ierr = MatSetValue(P, 1, 1, -1.0, INSERT_VALUES);CHKERRQ(ierr); 162*c4762a1bSJed Brown ierr = MatSetValue(P, 1, 2, 1.0, INSERT_VALUES);CHKERRQ(ierr); 163*c4762a1bSJed Brown 164*c4762a1bSJed Brown ierr = MatSetValue(P, 2, 0, 3.0, INSERT_VALUES);CHKERRQ(ierr); 165*c4762a1bSJed Brown ierr = MatSetValue(P, 2, 1, 0.0, INSERT_VALUES);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = MatSetValue(P, 2, 2, -3.0, INSERT_VALUES);CHKERRQ(ierr); 167*c4762a1bSJed Brown 168*c4762a1bSJed Brown ierr = MatAssemblyBegin(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 169*c4762a1bSJed Brown ierr = MatAssemblyEnd(P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 170*c4762a1bSJed Brown 171*c4762a1bSJed Brown /* compute C */ 172*c4762a1bSJed Brown ierr = MatPtAP(A, P, MAT_INITIAL_MATRIX, 1.0, &C);CHKERRQ(ierr); 173*c4762a1bSJed Brown 174*c4762a1bSJed Brown ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 175*c4762a1bSJed Brown ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 176*c4762a1bSJed Brown 177*c4762a1bSJed Brown /* compare results */ 178*c4762a1bSJed Brown /* 179*c4762a1bSJed Brown printf("C:\n"); 180*c4762a1bSJed Brown ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 181*c4762a1bSJed Brown 182*c4762a1bSJed Brown blitz::Array<double,2> actualC(cols, cols); 183*c4762a1bSJed Brown actualC = 0.0; 184*c4762a1bSJed Brown for (int i=0; i<cols; i++) { 185*c4762a1bSJed Brown for (int j=0; j<cols; j++) { 186*c4762a1bSJed Brown ierr = MatGetValues(C, 1, &i, 1, &j, &actualC(i,j)); 187*c4762a1bSJed Brown CHKERRQ(ierr); ; 188*c4762a1bSJed Brown } 189*c4762a1bSJed Brown } 190*c4762a1bSJed Brown blitz::Array<double,2> expectedC(cols, cols); 191*c4762a1bSJed Brown expectedC = 0.0; 192*c4762a1bSJed Brown 193*c4762a1bSJed Brown expectedC(0,0) = 10.0; 194*c4762a1bSJed Brown expectedC(0,1) = 2.0; 195*c4762a1bSJed Brown expectedC(0,2) = -9.0; 196*c4762a1bSJed Brown expectedC(0,3) = -1.0; 197*c4762a1bSJed Brown expectedC(1,0) = 2.0; 198*c4762a1bSJed Brown expectedC(1,1) = 5.0; 199*c4762a1bSJed Brown expectedC(1,2) = -1.0; 200*c4762a1bSJed Brown expectedC(1,3) = -2.0; 201*c4762a1bSJed Brown expectedC(2,0) = -9.0; 202*c4762a1bSJed Brown expectedC(2,1) = -1.0; 203*c4762a1bSJed Brown expectedC(2,2) = 10.0; 204*c4762a1bSJed Brown expectedC(2,3) = 0.0; 205*c4762a1bSJed Brown expectedC(3,0) = -1.0; 206*c4762a1bSJed Brown expectedC(3,1) = -2.0; 207*c4762a1bSJed Brown expectedC(3,2) = 0.0; 208*c4762a1bSJed Brown expectedC(3,3) = 1.0; 209*c4762a1bSJed Brown 210*c4762a1bSJed Brown int check = areBlitzArrays2NumericallyEqual(actualC,expectedC); 211*c4762a1bSJed Brown validateEqualsWithParams3(check, -1 , "testPTAPRectangular()", check, actualC(check), expectedC(check)); 212*c4762a1bSJed Brown */ 213*c4762a1bSJed Brown 214*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 215*c4762a1bSJed Brown ierr = MatDestroy(&P);CHKERRQ(ierr); 216*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 217*c4762a1bSJed Brown PetscFunctionReturn(0); 218*c4762a1bSJed Brown } 219*c4762a1bSJed Brown 220*c4762a1bSJed Brown /*TEST 221*c4762a1bSJed Brown 222*c4762a1bSJed Brown test: 223*c4762a1bSJed Brown 224*c4762a1bSJed Brown test: 225*c4762a1bSJed Brown suffix: 2 226*c4762a1bSJed Brown nsize: 2 227*c4762a1bSJed Brown args: -B_matmatmult_via nonscalable 228*c4762a1bSJed Brown 229*c4762a1bSJed Brown test: 230*c4762a1bSJed Brown suffix: 3 231*c4762a1bSJed Brown nsize: 2 232*c4762a1bSJed Brown output_file: output/ex93_2.out 233*c4762a1bSJed Brown 234*c4762a1bSJed Brown test: 235*c4762a1bSJed Brown suffix: 4 236*c4762a1bSJed Brown nsize: 2 237*c4762a1bSJed Brown args: -A_matptap_via scalable 238*c4762a1bSJed Brown output_file: output/ex93_2.out 239*c4762a1bSJed Brown 240*c4762a1bSJed Brown test: 241*c4762a1bSJed Brown suffix: btheap 242*c4762a1bSJed Brown args: -B_matmatmult_via btheap 243*c4762a1bSJed Brown output_file: output/ex93_1.out 244*c4762a1bSJed Brown 245*c4762a1bSJed Brown test: 246*c4762a1bSJed Brown suffix: heap 247*c4762a1bSJed Brown args: -B_matmatmult_via heap 248*c4762a1bSJed Brown output_file: output/ex93_1.out 249*c4762a1bSJed Brown 250*c4762a1bSJed Brown #HYPRE PtAP is broken for complex numbers 251*c4762a1bSJed Brown test: 252*c4762a1bSJed Brown suffix: hypre 253*c4762a1bSJed Brown nsize: 3 254*c4762a1bSJed Brown requires: hypre !complex 255*c4762a1bSJed Brown args: -B_matmatmult_via hypre -A_matptap_via hypre -test_hypre 256*c4762a1bSJed Brown 257*c4762a1bSJed Brown test: 258*c4762a1bSJed Brown suffix: llcondensed 259*c4762a1bSJed Brown args: -B_matmatmult_via llcondensed 260*c4762a1bSJed Brown output_file: output/ex93_1.out 261*c4762a1bSJed Brown 262*c4762a1bSJed Brown test: 263*c4762a1bSJed Brown suffix: scalable 264*c4762a1bSJed Brown args: -B_matmatmult_via scalable 265*c4762a1bSJed Brown output_file: output/ex93_1.out 266*c4762a1bSJed Brown 267*c4762a1bSJed Brown test: 268*c4762a1bSJed Brown suffix: scalable_fast 269*c4762a1bSJed Brown args: -B_matmatmult_via scalable_fast 270*c4762a1bSJed Brown output_file: output/ex93_1.out 271*c4762a1bSJed Brown 272*c4762a1bSJed Brown TEST*/ 273