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 PetscInt i,ia[2] = { 0, 2 }, ja[2] = { 0, 1 }, lda = 4; 11 PetscScalar a[2] = { 1.0, 1.0 }, *data; 12 PetscBool flg; 13 14 PetscCall(PetscInitialize(&argc, &argv, (char*)0,help)); 15 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 16 PetscCheck(size == 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must use 2 processors"); 17 18 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 19 PetscCall(MatSetType(A, MATMPIAIJ)); 20 PetscCall(MatSetSizes(A, 1, 1, 2, 2)); 21 PetscCall(MatMPIAIJSetPreallocationCSR(A, ia, ja, a)); 22 23 PetscCall(PetscCalloc1(4 * lda,&data)); 24 for (i = 0; i < 4; ++i) data[lda * i] = i * 1.0; 25 26 PetscCall(MatCreateDense(PETSC_COMM_WORLD, 1, PETSC_DECIDE, 2, 4, data, &B)); 27 PetscCall(MatSetOptionsPrefix(B,"b_")); 28 PetscCall(MatSetFromOptions(B)); 29 PetscCall(MatDenseSetLDA(B, lda)); 30 31 /* Test MatMatMult() */ 32 PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C)); 33 PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C)); 34 35 PetscCall(MatMatMultEqual(A,B,C,10,&flg)); 36 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C"); 37 38 /* Test user-provided mpidense matrix product */ 39 PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&C1)); 40 PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C1)); 41 PetscCall(MatMatMultEqual(A,B,C1,10,&flg)); 42 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C1"); 43 44 PetscCall(MatDestroy(&C1)); 45 PetscCall(MatDestroy(&C)); 46 47 /* Test MatTransposeMatMult() */ 48 PetscCall(MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C)); 49 PetscCall(MatTransposeMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C)); 50 51 PetscCall(MatTransposeMatMultEqual(A,B,C,10,&flg)); 52 PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatTransposeMatMult()"); 53 PetscCall(MatDestroy(&C)); 54 55 PetscCall(MatDestroy(&B)); 56 PetscCall(MatDestroy(&A)); 57 PetscCall(PetscFree(data)); 58 59 PetscCall(PetscFinalize()); 60 return 0; 61 } 62 63 /*TEST 64 65 test: 66 suffix: 1 67 nsize: 2 68 output_file: output/ex34.out 69 70 test: 71 suffix: 1_cuda 72 requires: cuda 73 nsize: 2 74 args: -b_mat_type mpidensecuda 75 output_file: output/ex34.out 76 TEST*/ 77