1c4762a1bSJed Brown static char help[] = "Test MatMatMult() and MatTransposeMatMult() for MPIAIJ and MPIDENSE matrices. \n\ 2ad16ce7aSStefano Zampini Sequential part of mpidense matrix allows changes made by MatDenseSetLDA(). \n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petsc.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc, char ** argv) 7c4762a1bSJed Brown { 8ad16ce7aSStefano Zampini Mat A, B, C, C1; 9c4762a1bSJed Brown PetscMPIInt size; 10c4762a1bSJed Brown PetscInt i,ia[2] = { 0, 2 }, ja[2] = { 0, 1 }, lda = 4; 11c4762a1bSJed Brown PetscScalar a[2] = { 1.0, 1.0 }, *data; 12c4762a1bSJed Brown PetscBool flg; 13c4762a1bSJed Brown 14*327415f7SBarry Smith PetscFunctionBeginUser; 159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char*)0,help)); 169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 17be096a46SBarry Smith PetscCheck(size == 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must use 2 processors"); 18c4762a1bSJed Brown 199566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 209566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATMPIAIJ)); 219566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, 1, 1, 2, 2)); 229566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocationCSR(A, ia, ja, a)); 23c4762a1bSJed Brown 249566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(4 * lda,&data)); 25c4762a1bSJed Brown for (i = 0; i < 4; ++i) data[lda * i] = i * 1.0; 26c4762a1bSJed Brown 279566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, 1, PETSC_DECIDE, 2, 4, data, &B)); 289566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(B,"b_")); 299566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 309566063dSJacob Faibussowitsch PetscCall(MatDenseSetLDA(B, lda)); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* Test MatMatMult() */ 339566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C)); 349566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C)); 35c4762a1bSJed Brown 369566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A,B,C,10,&flg)); 3728b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C"); 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Test user-provided mpidense matrix product */ 409566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&C1)); 419566063dSJacob Faibussowitsch PetscCall(MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C1)); 429566063dSJacob Faibussowitsch PetscCall(MatMatMultEqual(A,B,C1,10,&flg)); 4328b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C1"); 44c4762a1bSJed Brown 459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C1)); 469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* Test MatTransposeMatMult() */ 499566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C)); 509566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C)); 51c4762a1bSJed Brown 529566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMultEqual(A,B,C,10,&flg)); 5328b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatTransposeMatMult()"); 549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 589566063dSJacob Faibussowitsch PetscCall(PetscFree(data)); 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 61b122ec5aSJacob Faibussowitsch return 0; 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown /*TEST 65c4762a1bSJed Brown 66c4762a1bSJed Brown test: 67c4762a1bSJed Brown suffix: 1 68c4762a1bSJed Brown nsize: 2 69c4762a1bSJed Brown output_file: output/ex34.out 70a5225ed3SStefano Zampini 71a5225ed3SStefano Zampini test: 72a5225ed3SStefano Zampini suffix: 1_cuda 73a5225ed3SStefano Zampini requires: cuda 74a5225ed3SStefano Zampini nsize: 2 75a5225ed3SStefano Zampini args: -b_mat_type mpidensecuda 76a5225ed3SStefano Zampini output_file: output/ex34.out 77c4762a1bSJed Brown TEST*/ 78