Actual source code: ex78.c
2: 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\
3: Writes them using the PETSc sparse format.\n\
4: Note: I and J start at 1, not 0!\n\
5: Input parameters are:\n\
6: -Ain <filename> : input matrix in ascii format\n\
7: -bin <filename> : input rhs in ascii format\n\
8: -uin <filename> : input true solution in ascii format\n\
9: Run this program: ex33h -Ain Ain -bin bin -uin uin\n\n";
11: #include petscmat.h
15: int main(int argc,char **args)
16: {
17: Mat A;
18: Vec b,u,u_tmp;
19: char Ain[PETSC_MAX_PATH_LEN],bin[PETSC_MAX_PATH_LEN],uin[PETSC_MAX_PATH_LEN];
21: int m,n,nz,dummy,*col=0,*row=0; /* these are fscaned so kept as int */
22: PetscInt i,col_i,row_i,*nnz,*ib;
23: PetscScalar val_i,*work=0,mone=-1.0;
24: PetscReal res_norm,*val=0,*bval=0,*uval=0;
25: FILE *Afile,*bfile,*ufile;
26: PetscViewer view;
27: PetscTruth flg_A,flg_b,flg_u;
29: PetscInitialize(&argc,&args,(char *)0,help);
31: /* Read in matrix, rhs and exact solution from an ascii file */
32: PetscOptionsGetString(PETSC_NULL,"-Ain",Ain,PETSC_MAX_PATH_LEN-1,&flg_A);
33: if (flg_A){
34: PetscPrintf(PETSC_COMM_SELF,"\n Read matrix in ascii format ...\n");
35: PetscFOpen(PETSC_COMM_SELF,Ain,"r",&Afile);
36: fscanf(Afile,"%d %d %d\n",&m,&n,&nz);
37: printf("m: %d, n: %d, nz: %d \n", m,n,nz);
38: if (m != n) SETERRQ(PETSC_ERR_ARG_SIZ, "Number of rows, cols must be same for SBAIJ format\n");
40: VecCreateSeq(PETSC_COMM_SELF,n,&b);
41: VecCreateSeq(PETSC_COMM_SELF,n,&u);
43: PetscMalloc(m*sizeof(PetscInt),&nnz);
44: PetscMemzero(nnz,m*sizeof(PetscInt));
45: PetscMalloc(nz*sizeof(PetscReal),&val);
46: PetscMalloc(nz*sizeof(PetscInt),&row);
47: PetscMalloc(nz*sizeof(PetscInt),&col);
48: for (i=0; i<nz; i++) {
49: fscanf(Afile,"%d %d %le\n",row+i,col+i,(double*)val+i); /* modify according to data file! */
50: row[i]--; col[i]--; /* set index set starts at 0 */
51: nnz[row[i]]++;
52: }
53: printf("row[0]: %d, col[0]: %d, val: %g\n",row[0],col[0],val[0]);
54: printf("row[last]: %d, col: %d, val: %g\n",row[nz-1],col[nz-1],val[nz-1]);
55: fflush(stdout);
56: fclose(Afile);
57: MatCreateSeqBAIJ(PETSC_COMM_SELF,1,n,n,0,nnz,&A);
58: PetscFree(nnz);
59: }
61: PetscOptionsGetString(PETSC_NULL,"-bin",bin,PETSC_MAX_PATH_LEN-1,&flg_b);
62: if (flg_b){
63: PetscPrintf(PETSC_COMM_SELF,"\n Read rhs in ascii format ...\n");
64: PetscFOpen(PETSC_COMM_SELF,bin,"r",&bfile);
65: PetscMalloc(n*sizeof(PetscReal),&bval);
66: PetscMalloc(n*sizeof(PetscScalar),&work);
67: PetscMalloc(n*sizeof(PetscInt),&ib);
68: for (i=0; i<n; i++) {
69: /* fscanf(bfile,"%d %le\n",ib+i,bval+i); ib[i]--; */ /* modify according to data file! */
70: fscanf(bfile,"%d %le\n",&dummy,(double*)(bval+i)); ib[i] = i; /* modify according to data file! */
71: }
72: printf("bval[0]: %g, bval[%d]: %g\n",bval[0],(int)ib[n-1],bval[n-1]);
73: fflush(stdout);
74: fclose(bfile);
75: }
77: PetscOptionsGetString(PETSC_NULL,"-uin",uin,PETSC_MAX_PATH_LEN-1,&flg_u);
78: if (flg_u){
79: PetscPrintf(PETSC_COMM_SELF,"\n Read exact solution in ascii format ...\n");
80: PetscFOpen(PETSC_COMM_SELF,uin,"r",&ufile);
81: PetscMalloc(n*sizeof(PetscReal),&uval);
82: for (i=0; i<n; i++) {
83: fscanf(ufile,"%d %le\n",&dummy,(double*)(uval+i)); /* modify according to data file! */
84: }
85: printf("uval[0]: %g, uval[%d]: %g\n",uval[0], n-1, uval[n-1]);
86: fflush(stdout);
87: fclose(ufile);
88: }
90: if(flg_A){
91: /*
92: for (i=0; i<n; i++){
93: MatSetValues(A,1,&i,1,&i,&zero,INSERT_VALUES);
94: }
95: */
96: for (i=0; i<nz; i++) {
97: row_i =(int)row[i]; col_i =(int)col[i]; val_i = (PetscScalar)val[i];
98: MatSetValues(A,1,&row_i,1,&col_i,&val_i,INSERT_VALUES);
99: }
100: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
101: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
102: PetscFree(col);
103: PetscFree(val);
104: PetscFree(row);
105: }
106: if(flg_b){
107: for (i=0; i<n; i++) work[i]=(PetscScalar)bval[i];
108: VecSetValues(b,n,ib,work,INSERT_VALUES);
109: VecAssemblyBegin(b);
110: VecAssemblyEnd(b);
111: /* printf("b: \n"); VecView(b,PETSC_VIEWER_STDOUT_SELF); */
112: PetscFree(bval);
113: }
115: if(flg_u & flg_b){
116: for (i=0; i<n; i++) work[i]=(PetscScalar)uval[i];
117: VecSetValues(u,n,ib,work,INSERT_VALUES);
118: VecAssemblyBegin(u);
119: VecAssemblyEnd(u);
120: /* printf("u: \n"); VecView(u,PETSC_VIEWER_STDOUT_SELF); */
121: PetscFree(uval);
122: }
123:
124: if(flg_b) {
125: PetscFree(ib);
126: PetscFree(work);
127: }
128: /* Check accuracy of the data */
129: if (flg_A & flg_b & flg_u){
130: VecDuplicate(u,&u_tmp);
131: MatMult(A,u,u_tmp);
132: VecAXPY(&mone,b,u_tmp);
133: VecNorm(u_tmp,NORM_2,&res_norm);
134: printf("\n Accuracy of the reading data: | b - A*u |_2 : %g \n",res_norm);
136: /* Write the matrix, rhs and exact solution in Petsc binary file */
137: PetscPrintf(PETSC_COMM_SELF,"\n Write matrix and rhs in binary to 'matrix.dat' ...\n");
138: PetscViewerBinaryOpen(PETSC_COMM_SELF,"matrix.dat",PETSC_FILE_CREATE,&view);
139: MatView(A,view);
140: VecView(b,view);
141: VecView(u,view);
142: PetscViewerDestroy(view);
144: VecDestroy(b);
145: VecDestroy(u);
146: VecDestroy(u_tmp);
147:
148: MatDestroy(A);
149: }
151: PetscFinalize();
152: return 0;
153: }