Actual source code: ex34.c

petsc-3.13.6 2020-09-29
Report Typos and Errors
  1: static char help[] = "Test MatMatMult() and MatTransposeMatMult() for MPIAIJ and MPIDENSE matrices. \n\
  2:                       Sequential part of mpidense matrix allows changes made by MatSeqDenseSetLDA(). \n\n";

  4:  #include <petsc.h>

  6: int main(int argc, char ** argv)
  7: {
  8:   Mat            A, B, C, C1, seqB;
  9:   PetscMPIInt    size;
 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;

 15:   PetscInitialize(&argc, &argv, (char*)0,help);if (ierr) return ierr;
 16:   MPI_Comm_size(PETSC_COMM_WORLD, &size);
 17:   if (size != 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Must use 2 processors");

 19:   MatCreate(PETSC_COMM_WORLD, &A);
 20:   MatSetType(A, MATMPIAIJ);
 21:   MatSetSizes(A, 1, 1, 2, 2);
 22:   MatMPIAIJSetPreallocationCSR(A, ia, ja, a);

 24:   PetscCalloc1(4 * lda,&data);
 25:   for (i = 0; i < 4; ++i) data[lda * i] = i * 1.0;

 27:   MatCreateDense(PETSC_COMM_WORLD, 1, PETSC_DECIDE, 2, 4, data, &B);
 28:   MatDenseGetLocalMatrix(B, &seqB);
 29:   MatSeqDenseSetLDA(seqB, lda);

 31:   /* Test MatMatMult() */
 32:   MatMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C);
 33:   MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C);

 35:   MatMatMultEqual(A,B,C,10,&flg);
 36:   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C");

 38:   /* Test user-provided mpidense matrix product */
 39:   MatDuplicate(C,MAT_COPY_VALUES,&C1);
 40:   MatMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C1);
 41:   MatMatMultEqual(A,B,C1,10,&flg);
 42:   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatMatMult() for C1");

 44:   MatDestroy(&C1);
 45:   MatDestroy(&C);

 47:   /* Test MatTransposeMatMult() */
 48:   MatTransposeMatMult(A, B, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &C);
 49:   MatTransposeMatMult(A, B, MAT_REUSE_MATRIX, PETSC_DEFAULT, &C);

 51:   MatTransposeMatMultEqual(A,B,C,10,&flg);
 52:   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in MatTransposeMatMult()");
 53:   MatDestroy(&C);

 55:   MatDestroy(&B);
 56:   MatDestroy(&A);
 57:   PetscFree(data);

 59:   PetscFinalize();
 60:   return ierr;
 61: }

 63: /*TEST

 65:    test:
 66:       suffix: 1
 67:       nsize: 2
 68:       output_file: output/ex34.out
 69: TEST*/