Actual source code: ex47.c

  1: /*$Id: ex47.c,v 1.24 2001/08/07 21:30:08 bsmith Exp $*/

  3: static char help[] = "Tests the various routines in MatBAIJ format.n
  4: Input arguments are:n
  5:   -f <input_file> : file to load.  For a 5X5 example of the 5-pt. stencil,n
  6:                     use the file petsc/src/mat/examples/matbinary.exnn";

 8:  #include petscmat.h

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


 23:   PetscInitialize(&argc,&args,(char *)0,help);
 24: #if defined(PETSC_USE_COMPLEX)
 25:   SETERRQ(1,"This example does not work with complex numbers");
 26: #else
 27: 
 28:   PetscOptionsGetString(PETSC_NULL,"-f",file,127,PETSC_NULL);

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

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

 40:   /* Load the matrix as BAIJ format */
 41:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,PETSC_BINARY_RDONLY,&vc);
 42:   MatLoad(vc,MATSEQBAIJ,&C);
 43:   PetscViewerDestroy(vc);

 45:   MatGetSize(A,&m,&n);
 46:   MatGetSize(B,&m2,&n2);
 47:   if (m!=m2) SETERRQ(1,"Matrices are of different sixe. Cannot run this example");
 48: 
 49:   /* Test MatEqual() */
 50:   MatEqual(B,C,&tflg);
 51:   if (!tflg) SETERRQ(1,"MatEqual() failed");

 53:   /* Test MatGetDiagonal() */
 54:    VecCreateSeq(PETSC_COMM_SELF,m,&x);
 55:    VecCreateSeq(PETSC_COMM_SELF,m,&y);

 57:   MatGetDiagonal(A,x);
 58:   MatGetDiagonal(B,y);
 59: 
 60:   VecEqual(x,y,&tflg);
 61:   if (!tflg)  SETERRQ(1,"MatGetDiagonal() failed");

 63:   /* Test MatDiagonalScale() */
 64: 
 65:   PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&r);
 66:   VecSetRandom(r,x);
 67:   VecSetRandom(r,y);

 69:   MatDiagonalScale(A,x,y);
 70:   MatDiagonalScale(B,x,y);
 71:   MatMult(A,x,y);
 72:   VecNorm(y,NORM_2,&norm1);
 73:   MatMult(B,x,y);
 74:   VecNorm(y,NORM_2,&norm2);
 75:   rnorm = ((norm1-norm2)*100)/norm1;
 76:   if (rnorm<-0.1 || rnorm>0.01) {
 77:     PetscPrintf(PETSC_COMM_SELF,"Norm1=%e Norm2=%en",norm1,norm2);
 78:     SETERRQ(1,"MatDiagonalScale() failed");
 79:   }

 81:   /* Test MatGetRow()/ MatRestoreRow() */
 82:   for (ct=0; ct<100; ct++) {
 83:     PetscRandomGetValue(r,&rval);
 84:     row  = (int)(rval*m);
 85:     MatGetRow(A,row,&ncols1,&cols1,&vals1);
 86:     MatGetRow(B,row,&ncols2,&cols2,&vals2);
 87: 
 88:     for (i=0,j=0; i<ncols1 && j<ncols2; i++) {
 89:       while (cols2[j] != cols1[i]) j++;
 90:       if (vals1[i] != vals2[j]) SETERRQ(1,"MatGetRow() failed - vals incorrect.");
 91:     }
 92:     if (i<ncols1) SETERRQ(1,"MatGetRow() failed - cols incorrect");
 93: 
 94:     MatRestoreRow(A,row,&ncols1,&cols1,&vals1);
 95:     MatRestoreRow(B,row,&ncols2,&cols2,&vals2);
 96:   }
 97: 
 98:   MatDestroy(A);
 99:   MatDestroy(B);
100:   MatDestroy(C);
101:   VecDestroy(x);
102:   VecDestroy(y);
103:   PetscRandomDestroy(r);
104:   PetscFinalize();
105: #endif
106:   return 0;
107: }