Actual source code: ex78.c

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

  3: static char help[] = "Reads in a matrix in ASCII Matlab format (I,J,A), read in vectors rhs and exact_solu in ASCII format.n
  4: Writes them using the PETSc sparse format.n
  5: Note: I and J start at 1, not 0!n
  6: Input parameters are:n
  7:   -Ain  <filename> : input matrix in ascii formatn
  8:   -bin  <filename> : input rhs in ascii formatn
  9:   -uin  <filename> : input true solution in ascii formatn
 10: Run this program: ex33h -Ain Ain -bin bin -uin uinnn";

 12:  #include petscmat.h

 14: int main(int argc,char **args)
 15: {
 16:   Mat         A;
 17:   Vec         b,u,u_tmp;
 18:   char        Ain[128],bin[128],uin[128];
 19:   int         i,m,n,nz,ierr,*ib=0,col_i,row_i;
 20:   PetscScalar val_i,*work=0,mone=-1.0;
 21:   PetscReal   *col=0,*row=0,res_norm,*val=0,*bval=0,*uval=0;
 22:   FILE        *Afile,*bfile,*ufile;
 23:   PetscViewer view;
 24:   PetscTruth  flg_A,flg_b,flg_u;

 26:   PetscInitialize(&argc,&args,(char *)0,help);

 28:   /* Read in matrix, rhs and exact solution from an ascii file */
 29:   PetscOptionsGetString(PETSC_NULL,"-Ain",Ain,127,&flg_A);
 30:   if (flg_A){
 31:     PetscPrintf(PETSC_COMM_SELF,"n Read matrix in ascii format ...n");
 32:     PetscFOpen(PETSC_COMM_SELF,Ain,"r",&Afile);
 33:     fscanf(Afile,"%d %d %dn",&m,&n,&nz);
 34:     printf("m: %d, n: %d, nz: %d n", m,n,nz);
 35:     if (m != n) SETERRQ(PETSC_ERR_ARG_SIZ, "Number of rows, cols must be same for SBAIJ formatn");

 37:     MatCreateSeqBAIJ(PETSC_COMM_SELF,1,n,n,20,0,&A);
 38:     VecCreateSeq(PETSC_COMM_SELF,n,&b);
 39:     VecCreateSeq(PETSC_COMM_SELF,n,&u);

 41:     PetscMalloc(nz*sizeof(PetscReal),&val);
 42:     PetscMalloc(nz*sizeof(PetscReal),&row);
 43:     PetscMalloc(nz*sizeof(PetscReal),&col);
 44:     for (i=0; i<nz; i++) {
 45:       fscanf(Afile,"%le %le %len",row+i,col+i,val+i); /* modify according to data file! */
 46:       row[i]--; col[i]--;  /* set index set starts at 0 */
 47:     }
 48:     printf("row[0]: %g, col[0]: %g, val: %gn",row[0],col[0],val[0]);
 49:     printf("row[last]: %g, col: %g, val: %gn",row[nz-1],col[nz-1],val[nz-1]);
 50:     fclose(Afile);
 51:   }

 53:   PetscOptionsGetString(PETSC_NULL,"-bin",bin,127,&flg_b);
 54:   if (flg_b){
 55:     PetscPrintf(PETSC_COMM_SELF,"n Read rhs in ascii format ...n");
 56:     PetscFOpen(PETSC_COMM_SELF,bin,"r",&bfile);
 57:     PetscMalloc(n*sizeof(PetscReal),&bval);
 58:     PetscMalloc(n*sizeof(PetscScalar),&work);
 59:     PetscMalloc(n*sizeof(int),&ib);
 60:     for (i=0; i<n; i++) {
 61:       /* fscanf(bfile,"%d %len",ib+i,bval+i); ib[i]--;  */  /* modify according to data file! */
 62:       fscanf(bfile,"%len",bval+i); ib[i] = i;         /* modify according to data file! */
 63:     }
 64:     printf("bval[0]: %g, bval[%d]: %gn",bval[0],ib[n-1],bval[n-1]);
 65:     fclose(bfile);
 66:   }

 68:   PetscOptionsGetString(PETSC_NULL,"-uin",uin,127,&flg_u);
 69:   if (flg_u){
 70:     PetscPrintf(PETSC_COMM_SELF,"n Read exact solution in ascii format ...n");
 71:     PetscFOpen(PETSC_COMM_SELF,uin,"r",&ufile);
 72:     PetscMalloc(n*sizeof(PetscReal),&uval);
 73:     for (i=0; i<n; i++) {
 74:       fscanf(ufile,"  %len",uval+i);  /* modify according to data file! */
 75:     }
 76:     printf("uval[0]: %g, uval[%d]: %gn",uval[0], n-1, uval[n-1]);
 77:     fclose(ufile);
 78:   }

 80:   if(flg_A){
 81:     /*
 82:     for (i=0; i<n; i++){ 
 83:       MatSetValues(A,1,&i,1,&i,&zero,INSERT_VALUES); 
 84:     }
 85:     */
 86:     for (i=0; i<nz; i++) {
 87:       row_i =(int)row[i]; col_i =(int)col[i]; val_i = (PetscScalar)val[i];
 88:       MatSetValues(A,1,&row_i,1,&col_i,&val_i,INSERT_VALUES);
 89:     }
 90:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 91:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 92:     PetscFree(col);
 93:     PetscFree(val);
 94:     PetscFree(row);
 95:   }
 96:   if(flg_b){
 97:     for (i=0; i<n; i++) work[i]=(PetscScalar)bval[i];
 98:     VecSetValues(b,n,ib,work,INSERT_VALUES);
 99:     VecAssemblyBegin(b);
100:     VecAssemblyEnd(b);
101:     /* printf("b: n"); VecView(b,PETSC_VIEWER_STDOUT_SELF); */
102:     PetscFree(bval);
103:   }

105:   if(flg_u & flg_b){
106:     for (i=0; i<n; i++) work[i]=(PetscScalar)uval[i];
107:     VecSetValues(u,n,ib,work,INSERT_VALUES);
108:     VecAssemblyBegin(u);
109:     VecAssemblyEnd(u);
110:     /* printf("u: n"); VecView(u,PETSC_VIEWER_STDOUT_SELF); */
111:     PetscFree(uval);
112:   }
113: 
114:   if(flg_b) {
115:     PetscFree(ib);
116:     PetscFree(work);
117:   }
118:   /* Check accuracy of the data */
119:   if (flg_A & flg_b & flg_u){
120:     VecDuplicate(u,&u_tmp);
121:     MatMult(A,u,u_tmp);
122:     VecAXPY(&mone,b,u_tmp);
123:     VecNorm(u_tmp,NORM_2,&res_norm);
124:     printf("n Accuracy of the reading data: | b - A*u |_2 : %g n",res_norm);

126:   /* Write the matrix, rhs and exact solution in Petsc binary file */
127:     PetscPrintf(PETSC_COMM_SELF,"n Write matrix and rhs in binary to 'matrix.dat' ...n");
128:     PetscViewerBinaryOpen(PETSC_COMM_SELF,"matrix.dat",PETSC_BINARY_CREATE,&view);
129:     MatView(A,view);
130:     VecView(b,view);
131:     VecView(u,view);
132:     PetscViewerDestroy(view);

134:     VecDestroy(b);
135:     VecDestroy(u);
136:     VecDestroy(u_tmp);
137: 
138:     MatDestroy(A);
139:   }

141:   PetscFinalize();
142:   return 0;
143: }