1 static char help[] = "Test MatMatMult(), MatTranspose(), MatTransposeMatMult() for Dense and Elemental matrices.\n\n"; 2 /* 3 Example: 4 mpiexec -n <np> ./ex104 -mat_type elemental 5 */ 6 7 #include <petscmat.h> 8 9 int main(int argc,char **argv) 10 { 11 Mat A,B,C,D; 12 PetscInt i,M=10,N=5,j,nrows,ncols,am,an,rstart,rend; 13 PetscErrorCode ierr; 14 PetscRandom r; 15 PetscBool equal,Aiselemental; 16 PetscReal fill = 1.0; 17 IS isrows,iscols; 18 const PetscInt *rows,*cols; 19 PetscScalar *v,rval; 20 #if defined(PETSC_HAVE_ELEMENTAL) 21 PetscBool Test_MatMatMult=PETSC_TRUE; 22 #else 23 PetscBool Test_MatMatMult=PETSC_FALSE; 24 #endif 25 PetscMPIInt size; 26 27 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 28 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 29 30 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 31 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 32 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 33 CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); 34 CHKERRQ(MatSetType(A,MATDENSE)); 35 CHKERRQ(MatSetFromOptions(A)); 36 CHKERRQ(MatSetUp(A)); 37 CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&r)); 38 CHKERRQ(PetscRandomSetFromOptions(r)); 39 40 /* Set local matrix entries */ 41 CHKERRQ(MatGetOwnershipIS(A,&isrows,&iscols)); 42 CHKERRQ(ISGetLocalSize(isrows,&nrows)); 43 CHKERRQ(ISGetIndices(isrows,&rows)); 44 CHKERRQ(ISGetLocalSize(iscols,&ncols)); 45 CHKERRQ(ISGetIndices(iscols,&cols)); 46 CHKERRQ(PetscMalloc1(nrows*ncols,&v)); 47 for (i=0; i<nrows; i++) { 48 for (j=0; j<ncols; j++) { 49 CHKERRQ(PetscRandomGetValue(r,&rval)); 50 v[i*ncols+j] = rval; 51 } 52 } 53 CHKERRQ(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES)); 54 CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 55 CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 56 CHKERRQ(ISRestoreIndices(isrows,&rows)); 57 CHKERRQ(ISRestoreIndices(iscols,&cols)); 58 CHKERRQ(ISDestroy(&isrows)); 59 CHKERRQ(ISDestroy(&iscols)); 60 CHKERRQ(PetscRandomDestroy(&r)); 61 62 CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATELEMENTAL,&Aiselemental)); 63 64 /* Test MatCreateTranspose() and MatTranspose() */ 65 CHKERRQ(MatCreateTranspose(A,&C)); 66 CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&B)); /* B = A^T */ 67 CHKERRQ(MatMultEqual(C,B,10,&equal)); 68 PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"A^T*x != (x^T*A)^T"); 69 CHKERRQ(MatDestroy(&B)); 70 71 CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&B)); 72 if (!Aiselemental) { 73 CHKERRQ(MatTranspose(B,MAT_INPLACE_MATRIX,&B)); 74 CHKERRQ(MatMultEqual(C,B,10,&equal)); 75 PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"C*x != B*x"); 76 } 77 CHKERRQ(MatDestroy(&B)); 78 79 /* Test B = C*A for matrix type transpose and seqdense */ 80 if (size == 1 && !Aiselemental) { 81 CHKERRQ(MatMatMult(C,A,MAT_INITIAL_MATRIX,fill,&B)); 82 CHKERRQ(MatMatMultEqual(C,A,B,10,&equal)); 83 PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B != C*A for matrix type transpose and seqdense"); 84 CHKERRQ(MatDestroy(&B)); 85 } 86 CHKERRQ(MatDestroy(&C)); 87 88 /* Test MatMatMult() */ 89 if (Test_MatMatMult) { 90 #if !defined(PETSC_HAVE_ELEMENTAL) 91 PetscCheckFalse(size > 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"This test requires ELEMENTAL"); 92 #endif 93 CHKERRQ(MatTranspose(A,MAT_INITIAL_MATRIX,&B)); /* B = A^T */ 94 CHKERRQ(MatMatMult(B,A,MAT_INITIAL_MATRIX,fill,&C)); /* C = B*A = A^T*A */ 95 CHKERRQ(MatMatMult(B,A,MAT_REUSE_MATRIX,fill,&C)); 96 CHKERRQ(MatMatMultEqual(B,A,C,10,&equal)); 97 PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"B*A*x != C*x"); 98 99 /* Test MatDuplicate for matrix product */ 100 CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&D)); 101 102 CHKERRQ(MatDestroy(&D)); 103 CHKERRQ(MatDestroy(&C)); 104 CHKERRQ(MatDestroy(&B)); 105 } 106 107 /* Test MatTransposeMatMult() */ 108 if (!Aiselemental) { 109 CHKERRQ(MatTransposeMatMult(A,A,MAT_INITIAL_MATRIX,fill,&D)); /* D = A^T*A */ 110 CHKERRQ(MatTransposeMatMult(A,A,MAT_REUSE_MATRIX,fill,&D)); 111 CHKERRQ(MatTransposeMatMultEqual(A,A,D,10,&equal)); 112 PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x"); 113 114 /* Test MatDuplicate for matrix product */ 115 CHKERRQ(MatDuplicate(D,MAT_COPY_VALUES,&C)); 116 CHKERRQ(MatDestroy(&C)); 117 CHKERRQ(MatDestroy(&D)); 118 119 /* Test D*x = A^T*C*A*x, where C is in AIJ format */ 120 CHKERRQ(MatGetLocalSize(A,&am,&an)); 121 CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 122 if (size == 1) { 123 CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,am,am)); 124 } else { 125 CHKERRQ(MatSetSizes(C,am,am,PETSC_DECIDE,PETSC_DECIDE)); 126 } 127 CHKERRQ(MatSetFromOptions(C)); 128 CHKERRQ(MatSetUp(C)); 129 CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend)); 130 v[0] = 1.0; 131 for (i=rstart; i<rend; i++) { 132 CHKERRQ(MatSetValues(C,1,&i,1,&i,v,INSERT_VALUES)); 133 } 134 CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 135 CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 136 137 /* B = C*A, D = A^T*B */ 138 CHKERRQ(MatMatMult(C,A,MAT_INITIAL_MATRIX,1.0,&B)); 139 CHKERRQ(MatTransposeMatMult(A,B,MAT_INITIAL_MATRIX,fill,&D)); 140 CHKERRQ(MatTransposeMatMultEqual(A,B,D,10,&equal)); 141 PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*B*x"); 142 143 CHKERRQ(MatDestroy(&D)); 144 CHKERRQ(MatDestroy(&C)); 145 CHKERRQ(MatDestroy(&B)); 146 } 147 148 /* Test MatMatTransposeMult() */ 149 if (!Aiselemental) { 150 PetscReal diff, scale; 151 PetscInt am, an, aM, aN; 152 153 CHKERRQ(MatGetLocalSize(A, &am, &an)); 154 CHKERRQ(MatGetSize(A, &aM, &aN)); 155 CHKERRQ(MatCreateDense(PetscObjectComm((PetscObject)A),PETSC_DECIDE, an, aM + 10, aN, NULL, &B)); 156 CHKERRQ(MatSetRandom(B, NULL)); 157 CHKERRQ(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 158 CHKERRQ(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 159 CHKERRQ(MatMatTransposeMult(A,B,MAT_INITIAL_MATRIX,fill,&D)); /* D = A*A^T */ 160 161 /* Test MatDuplicate for matrix product */ 162 CHKERRQ(MatDuplicate(D,MAT_COPY_VALUES,&C)); 163 164 CHKERRQ(MatMatTransposeMult(A,B,MAT_REUSE_MATRIX,fill,&D)); 165 CHKERRQ(MatAXPY(C, -1., D, SAME_NONZERO_PATTERN)); 166 167 CHKERRQ(MatNorm(C, NORM_FROBENIUS, &diff)); 168 CHKERRQ(MatNorm(D, NORM_FROBENIUS, &scale)); 169 PetscCheckFalse(diff > PETSC_SMALL * scale,PetscObjectComm((PetscObject)D), PETSC_ERR_PLIB, "MatMatTransposeMult() differs between MAT_INITIAL_MATRIX and MAT_REUSE_MATRIX"); 170 CHKERRQ(MatDestroy(&C)); 171 172 CHKERRQ(MatMatTransposeMultEqual(A,B,D,10,&equal)); 173 PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_PLIB,"D*x != A^T*A*x"); 174 CHKERRQ(MatDestroy(&D)); 175 CHKERRQ(MatDestroy(&B)); 176 177 } 178 179 CHKERRQ(MatDestroy(&A)); 180 CHKERRQ(PetscFree(v)); 181 ierr = PetscFinalize(); 182 return ierr; 183 } 184 185 /*TEST 186 187 test: 188 output_file: output/ex104.out 189 190 test: 191 suffix: 2 192 nsize: 2 193 output_file: output/ex104.out 194 195 test: 196 suffix: 3 197 nsize: 4 198 output_file: output/ex104.out 199 args: -M 23 -N 31 200 201 test: 202 suffix: 4 203 nsize: 4 204 output_file: output/ex104.out 205 args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via cyclic 206 207 test: 208 suffix: 5 209 nsize: 4 210 output_file: output/ex104.out 211 args: -M 23 -N 31 -matmattransmult_mpidense_mpidense_via allgatherv 212 213 test: 214 suffix: 6 215 args: -mat_type elemental 216 requires: elemental 217 output_file: output/ex104.out 218 219 test: 220 suffix: 7 221 nsize: 2 222 args: -mat_type elemental 223 requires: elemental 224 output_file: output/ex104.out 225 226 TEST*/ 227