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.ex\n\n";

 8:  #include petscmat.h

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


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

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

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

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

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

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

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

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

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

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