Actual source code: ex100.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: static char help[] = "Tests vatious routines in MatMAIJ format.\n";

  4: #include <petscmat.h>
  5: #define IMAX 15
  8: int main(int argc,char **args)
  9: {
 10:   Mat            A,B,MA;
 11:   PetscViewer    fd;
 12:   char           file[PETSC_MAX_PATH_LEN];
 13:   PetscInt       m,n,M,N,dof=1;
 14:   PetscMPIInt    rank,size;
 16:   PetscBool      flg;

 18:   PetscInitialize(&argc,&args,(char*)0,help);
 19:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 20:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 22: #if defined(PETSC_USE_COMPLEX)
 23:   SETERRQ(PETSC_COMM_WORLD,1,"This example does not work with complex numbers");
 24: #else

 26:   /* Load aij matrix A */
 27:   PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg);
 28:   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
 29:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);
 30:   MatCreate(PETSC_COMM_WORLD,&A);
 31:   MatLoad(A,fd);
 32:   PetscViewerDestroy(&fd);

 34:   /* Get dof, then create maij matrix MA */
 35:   PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL);
 36:   MatCreateMAIJ(A,dof,&MA);
 37:   MatGetLocalSize(MA,&m,&n);
 38:   MatGetSize(MA,&M,&N);

 40:   if (size == 1) {
 41:     MatConvert(MA,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
 42:   } else {
 43:     MatConvert(MA,MATMPIAIJ,MAT_INITIAL_MATRIX,&B);
 44:   }

 46:   /* Test MatMult() */
 47:   MatMultEqual(MA,B,10,&flg);
 48:   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error: MatMul() for MAIJ matrix");
 49:   /* Test MatMultAdd() */
 50:   MatMultAddEqual(MA,B,10,&flg);
 51:   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error: MatMulAdd() for MAIJ matrix");

 53:   /* Test MatMultTranspose() */
 54:   MatMultTransposeEqual(MA,B,10,&flg);
 55:   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error: MatMulAdd() for MAIJ matrix");

 57:   /* Test MatMultTransposeAdd() */
 58:   MatMultTransposeAddEqual(MA,B,10,&flg);
 59:   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Error: MatMulTransposeAdd() for MAIJ matrix");

 61:   MatDestroy(&MA);
 62:   MatDestroy(&A);
 63:   MatDestroy(&B);
 64:   PetscFinalize();
 65: #endif
 66:   return 0;
 67: }