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 PetscErrorCode ierr; 11c4762a1bSJed Brown PetscInt i,ia[2] = { 0, 2 }, ja[2] = { 0, 1 }, lda = 4; 12c4762a1bSJed Brown PetscScalar a[2] = { 1.0, 1.0 }, *data; 13c4762a1bSJed Brown PetscBool flg; 14c4762a1bSJed Brown 15c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, (char*)0,help);if (ierr) return ierr; 16ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr); 17*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 2,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Must use 2 processors"); 18c4762a1bSJed Brown 19c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD, &A);CHKERRQ(ierr); 20c4762a1bSJed Brown ierr = MatSetType(A, MATMPIAIJ);CHKERRQ(ierr); 21c4762a1bSJed Brown ierr = MatSetSizes(A, 1, 1, 2, 2);CHKERRQ(ierr); 22c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocationCSR(A, ia, ja, a);CHKERRQ(ierr); 23c4762a1bSJed Brown 24c4762a1bSJed Brown ierr = PetscCalloc1(4 * lda,&data);CHKERRQ(ierr); 25c4762a1bSJed Brown for (i = 0; i < 4; ++i) data[lda * i] = i * 1.0; 26c4762a1bSJed Brown 27c4762a1bSJed Brown ierr = MatCreateDense(PETSC_COMM_WORLD, 1, PETSC_DECIDE, 2, 4, data, &B);CHKERRQ(ierr); 28a5225ed3SStefano Zampini ierr = MatSetOptionsPrefix(B,"b_");CHKERRQ(ierr); 29a5225ed3SStefano Zampini ierr = MatSetFromOptions(B);CHKERRQ(ierr); 30ad16ce7aSStefano Zampini ierr = MatDenseSetLDA(B, lda);CHKERRQ(ierr); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* Test MatMatMult() */ 33c4762a1bSJed Brown ierr = MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C);CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C);CHKERRQ(ierr); 35c4762a1bSJed Brown 36c4762a1bSJed Brown ierr = MatMatMultEqual(A,B,C,10,&flg);CHKERRQ(ierr); 37*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C"); 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Test user-provided mpidense matrix product */ 40c4762a1bSJed Brown ierr = MatDuplicate(C,MAT_COPY_VALUES,&C1);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C1);CHKERRQ(ierr); 42c4762a1bSJed Brown ierr = MatMatMultEqual(A,B,C1,10,&flg);CHKERRQ(ierr); 43*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C1"); 44c4762a1bSJed Brown 45c4762a1bSJed Brown ierr = MatDestroy(&C1);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* Test MatTransposeMatMult() */ 49c4762a1bSJed Brown ierr = MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = MatTransposeMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C);CHKERRQ(ierr); 51c4762a1bSJed Brown 52c4762a1bSJed Brown ierr = MatTransposeMatMultEqual(A,B,C,10,&flg);CHKERRQ(ierr); 53*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatTransposeMatMult()"); 54c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 55c4762a1bSJed Brown 56c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 57c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = PetscFree(data);CHKERRQ(ierr); 59c4762a1bSJed Brown 60c4762a1bSJed Brown ierr = PetscFinalize(); 61c4762a1bSJed Brown return ierr; 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