Actual source code: ex47.c

petsc-3.8.4 2018-03-24
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: #if defined(PETSC_USE_COMPLEX)
 26:   SETERRQ(PETSC_COMM_WORLD,1,"This example does not work with complex numbers");
 27: #else
 28:   PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,NULL);

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

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

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

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

 55:   /* Test MatEqual() */
 56:   MatEqual(B,C,&tflg);
 57:   if (!tflg) SETERRQ(PETSC_COMM_SELF,1,"MatEqual() failed");

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

 63:   MatGetDiagonal(A,x);
 64:   MatGetDiagonal(B,y);

 66:   VecEqual(x,y,&tflg);
 67:   if (!tflg) SETERRQ(PETSC_COMM_SELF,1,"MatGetDiagonal() failed");

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

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

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

 93:     for (i=0,j=0; i<ncols1 && j<ncols2; i++) {
 94:       while (cols2[j] != cols1[i]) j++;
 95:       if (vals1[i] != vals2[j]) SETERRQ(PETSC_COMM_SELF,1,"MatGetRow() failed - vals incorrect.");
 96:     }
 97:     if (i<ncols1) SETERRQ(PETSC_COMM_SELF,1,"MatGetRow() failed - cols incorrect");

 99:     MatRestoreRow(A,row,&ncols1,&cols1,&vals1);
100:     MatRestoreRow(B,row,&ncols2,&cols2,&vals2);
101:   }

103:   MatDestroy(&A);
104:   MatDestroy(&B);
105:   MatDestroy(&C);
106:   VecDestroy(&x);
107:   VecDestroy(&y);
108:   PetscRandomDestroy(&r);
109: #endif
110:   PetscFinalize();
111:   return ierr;
112: }