1 static char help[] = "Test MatMatMult() and MatTransposeMatMult() for MPIAIJ and MPIDENSE matrices. \n\ 2 Sequential part of mpidense matrix allows changes made by MatDenseSetLDA(). \n\n"; 3 4 #include <petsc.h> 5 6 int main(int argc, char ** argv) 7 { 8 Mat A, B, C, C1; 9 PetscMPIInt size; 10 PetscErrorCode ierr; 11 PetscInt i,ia[2] = { 0, 2 }, ja[2] = { 0, 1 }, lda = 4; 12 PetscScalar a[2] = { 1.0, 1.0 }, *data; 13 PetscBool flg; 14 15 ierr = PetscInitialize(&argc, &argv, (char*)0,help);if (ierr) return ierr; 16 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 17 PetscCheckFalse(size != 2,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Must use 2 processors"); 18 19 CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A)); 20 CHKERRQ(MatSetType(A, MATMPIAIJ)); 21 CHKERRQ(MatSetSizes(A, 1, 1, 2, 2)); 22 CHKERRQ(MatMPIAIJSetPreallocationCSR(A, ia, ja, a)); 23 24 CHKERRQ(PetscCalloc1(4 * lda,&data)); 25 for (i = 0; i < 4; ++i) data[lda * i] = i * 1.0; 26 27 CHKERRQ(MatCreateDense(PETSC_COMM_WORLD, 1, PETSC_DECIDE, 2, 4, data, &B)); 28 CHKERRQ(MatSetOptionsPrefix(B,"b_")); 29 CHKERRQ(MatSetFromOptions(B)); 30 CHKERRQ(MatDenseSetLDA(B, lda)); 31 32 /* Test MatMatMult() */ 33 CHKERRQ(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C)); 34 CHKERRQ(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C)); 35 36 CHKERRQ(MatMatMultEqual(A,B,C,10,&flg)); 37 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C"); 38 39 /* Test user-provided mpidense matrix product */ 40 CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&C1)); 41 CHKERRQ(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C1)); 42 CHKERRQ(MatMatMultEqual(A,B,C1,10,&flg)); 43 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C1"); 44 45 CHKERRQ(MatDestroy(&C1)); 46 CHKERRQ(MatDestroy(&C)); 47 48 /* Test MatTransposeMatMult() */ 49 CHKERRQ(MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C)); 50 CHKERRQ(MatTransposeMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C)); 51 52 CHKERRQ(MatTransposeMatMultEqual(A,B,C,10,&flg)); 53 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatTransposeMatMult()"); 54 CHKERRQ(MatDestroy(&C)); 55 56 CHKERRQ(MatDestroy(&B)); 57 CHKERRQ(MatDestroy(&A)); 58 CHKERRQ(PetscFree(data)); 59 60 ierr = PetscFinalize(); 61 return ierr; 62 } 63 64 /*TEST 65 66 test: 67 suffix: 1 68 nsize: 2 69 output_file: output/ex34.out 70 71 test: 72 suffix: 1_cuda 73 requires: cuda 74 nsize: 2 75 args: -b_mat_type mpidensecuda 76 output_file: output/ex34.out 77 TEST*/ 78