Actual source code: ex47.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2: static char help[] = "Tests the various routines in MatBAIJ format.\n\
  3: Input arguments are:\n\
  4:   -f <input_file> : file to load.  For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\n";

  6:  #include <petscmat.h>

  8: int main(int argc,char **args)
  9: {
 10:   Mat               A,B,C;
 11:   PetscViewer       va,vb,vc;
 12:   Vec               x,y;
 13:   PetscErrorCode    ierr;
 14:   PetscInt          i,j,row,m,n,ncols1,ncols2,ct,m2,n2;
 15:   const PetscInt    *cols1,*cols2;
 16:   char              file[PETSC_MAX_PATH_LEN];
 17:   PetscBool         tflg;
 18:   PetscScalar       rval;
 19:   const PetscScalar *vals1,*vals2;
 20:   PetscReal         norm1,norm2,rnorm;
 21:   PetscRandom       r;


 24:   PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
 25:   PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,NULL);

 27:   /* Load the matrix as AIJ format */
 28:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&va);
 29:   MatCreate(PETSC_COMM_WORLD,&A);
 30:   MatSetType(A,MATSEQAIJ);
 31:   MatLoad(A,va);
 32:   PetscViewerDestroy(&va);

 34:   /* Load the matrix as BAIJ format */
 35:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vb);
 36:   MatCreate(PETSC_COMM_WORLD,&B);
 37:   MatSetType(B,MATSEQBAIJ);
 38:   MatLoad(B,vb);
 39:   PetscViewerDestroy(&vb);

 41:   /* Load the matrix as BAIJ format */
 42:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vc);
 43:   MatCreate(PETSC_COMM_WORLD,&C);
 44:   MatSetType(C,MATSEQBAIJ);
 45:   MatSetFromOptions(C);
 46:   MatLoad(C,vc);
 47:   PetscViewerDestroy(&vc);

 49:   MatGetSize(A,&m,&n);
 50:   MatGetSize(B,&m2,&n2);
 51:   if (m!=m2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices are of different size. Cannot run this example");

 53:   /* Test MatEqual() */
 54:   MatEqual(B,C,&tflg);
 55:   if (!tflg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatEqual() failed");

 57:   /* Test MatGetDiagonal() */
 58:   VecCreateSeq(PETSC_COMM_SELF,m,&x);
 59:   VecCreateSeq(PETSC_COMM_SELF,m,&y);

 61:   MatGetDiagonal(A,x);
 62:   MatGetDiagonal(B,y);

 64:   VecEqual(x,y,&tflg);
 65:   if (!tflg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetDiagonal() failed");

 67:   /* Test MatDiagonalScale() */
 68:   PetscRandomCreate(PETSC_COMM_SELF,&r);
 69:   PetscRandomSetFromOptions(r);
 70:   VecSetRandom(x,r);
 71:   VecSetRandom(y,r);

 73:   MatDiagonalScale(A,x,y);
 74:   MatDiagonalScale(B,x,y);
 75:   MatMult(A,x,y);
 76:   VecNorm(y,NORM_2,&norm1);
 77:   MatMult(B,x,y);
 78:   VecNorm(y,NORM_2,&norm2);
 79:   rnorm = ((norm1-norm2)*100)/norm1;
 80:   if (rnorm<-0.1 || rnorm>0.01) SETERRQ2(PETSC_COMM_SELF,1,"MatDiagonalScale() failed Norm1 %g Norm2 %g",(double)norm1,(double)norm2);

 82:   /* Test MatGetRow()/ MatRestoreRow() */
 83:   for (ct=0; ct<100; ct++) {
 84:     PetscRandomGetValue(r,&rval);
 85:     row  = (int)(rval*m);
 86:     MatGetRow(A,row,&ncols1,&cols1,&vals1);
 87:     MatGetRow(B,row,&ncols2,&cols2,&vals2);

 89:     for (i=0,j=0; i<ncols1 && j<ncols2; i++) {
 90:       while (cols2[j] != cols1[i]) j++;
 91:       if (vals1[i] != vals2[j]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - vals incorrect.");
 92:     }
 93:     if (i<ncols1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - cols incorrect");

 95:     MatRestoreRow(A,row,&ncols1,&cols1,&vals1);
 96:     MatRestoreRow(B,row,&ncols2,&cols2,&vals2);
 97:   }

 99:   MatDestroy(&A);
100:   MatDestroy(&B);
101:   MatDestroy(&C);
102:   VecDestroy(&x);
103:   VecDestroy(&y);
104:   PetscRandomDestroy(&r);
105:   PetscFinalize();
106:   return ierr;
107: }

109: /*TEST

111:    build:
112:       requires: !complex

114:    test:
115:       args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -mat_block_size 5
116:       requires: !complex double datafilespath !define(PETSC_USE_64BIT_INDICES)

118: TEST*/