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