Actual source code: ex100.c

  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;
 15:   PetscErrorCode    ierr;
 16:   PetscTruth        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(1,"This example does not work with complex numbers");
 24: #else

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

 33:   /* Get dof, then create maij matrix MA */
 34:   PetscOptionsGetInt(PETSC_NULL,"-dof",&dof,PETSC_NULL);
 35:   MatCreateMAIJ(A,dof,&MA);
 36:   MatGetLocalSize(MA,&m,&n);
 37:   MatGetSize(MA,&M,&N);
 38: 
 39:   if (size == 1){
 40:     MatConvert(MA,MATSEQAIJ,MAT_INITIAL_MATRIX,&B);
 41:   } else {
 42:     MatConvert(MA,MATMPIAIJ,MAT_INITIAL_MATRIX,&B);
 43:   }

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

 56:   /* Test MatMultTranspose() */
 57:   MatMultTransposeEqual(MA,B,10,&flg);
 58:   if (!flg){
 59:     SETERRQ(PETSC_ERR_CONV_FAILED,"Error: MatMulAdd() for MAIJ matrix");
 60:   }
 61: 
 62:   /* Test MatMultTransposeAdd() */
 63:    MatMultTransposeAddEqual(MA,B,10,&flg);
 64:   if (!flg){
 65:     SETERRQ(PETSC_ERR_CONV_FAILED,"Error: MatMulTransposeAdd() for MAIJ matrix");
 66:   }

 68:   MatDestroy(MA);
 69:   MatDestroy(A);
 70:   MatDestroy(B);
 71:   PetscFinalize();
 72: #endif
 73:   return 0;
 74: }