1 static char help[] = "Test MatMatMult() for AIJ and Dense matrices.\n\n"; 2 3 #include <petscmat.h> 4 5 int main(int argc,char **argv) 6 { 7 Mat A,B,C,D,AT; 8 PetscInt i,M,N,Istart,Iend,n=7,j,J,Ii,m=8,am,an; 9 PetscScalar v; 10 PetscRandom r; 11 PetscBool equal=PETSC_FALSE,flg; 12 PetscReal fill = 1.0,norm; 13 PetscMPIInt size; 14 15 CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 16 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 17 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 18 CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-fill",&fill,NULL)); 19 20 CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&r)); 21 CHKERRQ(PetscRandomSetFromOptions(r)); 22 23 /* Create a aij matrix A */ 24 M = N = m*n; 25 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 26 CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); 27 CHKERRQ(MatSetType(A,MATAIJ)); 28 CHKERRQ(MatSetFromOptions(A)); 29 CHKERRQ(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL)); 30 CHKERRQ(MatSeqAIJSetPreallocation(A,5,NULL)); 31 32 CHKERRQ(MatGetOwnershipRange(A,&Istart,&Iend)); 33 am = Iend - Istart; 34 for (Ii=Istart; Ii<Iend; Ii++) { 35 v = -1.0; i = Ii/n; j = Ii - i*n; 36 if (i>0) {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 37 if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 38 if (j>0) {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 39 if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 40 v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 41 } 42 CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 43 CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 44 45 /* Create a dense matrix B */ 46 CHKERRQ(MatGetLocalSize(A,&am,&an)); 47 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 48 CHKERRQ(MatSetSizes(B,an,PETSC_DECIDE,PETSC_DECIDE,M)); 49 CHKERRQ(MatSetType(B,MATDENSE)); 50 CHKERRQ(MatSeqDenseSetPreallocation(B,NULL)); 51 CHKERRQ(MatMPIDenseSetPreallocation(B,NULL)); 52 CHKERRQ(MatSetFromOptions(B)); 53 CHKERRQ(MatSetRandom(B,r)); 54 CHKERRQ(PetscRandomDestroy(&r)); 55 CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 56 CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 57 58 /* Test reuse of user-provided dense C (unassembled) -- not recommended usage */ 59 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 60 CHKERRQ(MatSetType(C,MATDENSE)); 61 CHKERRQ(MatSetSizes(C,am,PETSC_DECIDE,PETSC_DECIDE,N)); 62 CHKERRQ(MatSetFromOptions(C)); 63 CHKERRQ(MatSetUp(C)); 64 CHKERRQ(MatZeroEntries(C)); 65 CHKERRQ(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C)); 66 CHKERRQ(MatNorm(C,NORM_INFINITY,&norm)); 67 CHKERRQ(MatDestroy(&C)); 68 69 /* Test C = A*B (aij*dense) */ 70 CHKERRQ(MatMatMult(A,B,MAT_INITIAL_MATRIX,fill,&C)); 71 CHKERRQ(MatMatMult(A,B,MAT_REUSE_MATRIX,fill,&C)); 72 73 /* Test developer API */ 74 CHKERRQ(MatProductCreate(A,B,NULL,&D)); 75 CHKERRQ(MatProductSetType(D,MATPRODUCT_AB)); 76 CHKERRQ(MatProductSetAlgorithm(D,"default")); 77 CHKERRQ(MatProductSetFill(D,fill)); 78 CHKERRQ(MatProductSetFromOptions(D)); 79 CHKERRQ(MatProductSymbolic(D)); 80 for (i=0; i<2; i++) { 81 CHKERRQ(MatProductNumeric(D)); 82 } 83 CHKERRQ(MatEqual(C,D,&equal)); 84 PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"C != D"); 85 CHKERRQ(MatDestroy(&D)); 86 87 /* Test D = AT*B (transpose(aij)*dense) */ 88 CHKERRQ(MatCreateTranspose(A,&AT)); 89 CHKERRQ(MatMatMult(AT,B,MAT_INITIAL_MATRIX,fill,&D)); 90 CHKERRQ(MatMatMultEqual(AT,B,D,10,&equal)); 91 PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != AT*B (transpose(aij)*dense)"); 92 CHKERRQ(MatDestroy(&D)); 93 CHKERRQ(MatDestroy(&AT)); 94 95 /* Test D = C*A (dense*aij) */ 96 CHKERRQ(MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&D)); 97 CHKERRQ(MatMatMult(C,A,MAT_REUSE_MATRIX,fill,&D)); 98 CHKERRQ(MatMatMultEqual(C,A,D,10,&equal)); 99 PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != C*A (dense*aij)"); 100 CHKERRQ(MatDestroy(&D)); 101 102 /* Test D = A*C (aij*dense) */ 103 CHKERRQ(MatMatMult(A,C,MAT_INITIAL_MATRIX,fill,&D)); 104 CHKERRQ(MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&D)); 105 CHKERRQ(MatMatMultEqual(A,C,D,10,&equal)); 106 PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != A*C (aij*dense)"); 107 CHKERRQ(MatDestroy(&D)); 108 109 /* Test D = B*C (dense*dense) */ 110 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 111 if (size == 1) { 112 CHKERRQ(MatMatMult(B,C,MAT_INITIAL_MATRIX,fill,&D)); 113 CHKERRQ(MatMatMult(B,C,MAT_REUSE_MATRIX,fill,&D)); 114 CHKERRQ(MatMatMultEqual(B,C,D,10,&equal)); 115 PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C (dense*dense)"); 116 CHKERRQ(MatDestroy(&D)); 117 } 118 119 /* Test D = B*C^T (dense*dense) */ 120 CHKERRQ(MatMatTransposeMult(B,C,MAT_INITIAL_MATRIX,fill,&D)); 121 CHKERRQ(MatMatTransposeMult(B,C,MAT_REUSE_MATRIX,fill,&D)); 122 CHKERRQ(MatMatTransposeMultEqual(B,C,D,10,&equal)); 123 PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"D != B*C^T (dense*dense)"); 124 CHKERRQ(MatDestroy(&D)); 125 126 /* Test MatProductCreateWithMat() and reuse C and B for B = A*C */ 127 flg = PETSC_FALSE; 128 CHKERRQ(PetscOptionsHasName(NULL,NULL,"-test_userAPI",&flg)); 129 if (flg) { 130 /* user driver */ 131 CHKERRQ(MatMatMult(A,C,MAT_REUSE_MATRIX,fill,&B)); 132 } else { 133 /* clear internal data structures related with previous products to avoid circular references */ 134 CHKERRQ(MatProductClear(A)); 135 CHKERRQ(MatProductClear(B)); 136 CHKERRQ(MatProductClear(C)); 137 CHKERRQ(MatProductCreateWithMat(A,C,NULL,B)); 138 CHKERRQ(MatProductSetType(B,MATPRODUCT_AB)); 139 CHKERRQ(MatProductSetFromOptions(B)); 140 CHKERRQ(MatProductSymbolic(B)); 141 CHKERRQ(MatProductNumeric(B)); 142 } 143 144 CHKERRQ(MatDestroy(&C)); 145 CHKERRQ(MatDestroy(&B)); 146 CHKERRQ(MatDestroy(&A)); 147 CHKERRQ(PetscFinalize()); 148 return 0; 149 } 150 151 /*TEST 152 153 test: 154 args: -M 10 -N 10 155 output_file: output/ex109.out 156 157 test: 158 suffix: 2 159 nsize: 3 160 output_file: output/ex109.out 161 162 test: 163 suffix: 3 164 nsize: 2 165 args: -matmattransmult_mpidense_mpidense_via cyclic 166 output_file: output/ex109.out 167 168 test: 169 suffix: 4 170 args: -test_userAPI 171 output_file: output/ex109.out 172 173 test: 174 suffix: 5 175 nsize: 3 176 args: -test_userAPI 177 output_file: output/ex109.out 178 179 TEST*/ 180