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