Actual source code: ex78.c
petsc-3.8.4 2018-03-24
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, use -noshift if indices in file start with zero!\n\
5: Input parameters are:\n\
6: -Ain <filename> : input matrix in ascii format\n\
7: -rhs <filename> : input rhs in ascii format\n\
8: -solu <filename> : input true solution in ascii format\n\\n";
10: /*
11: Example: ./ex78 -Ain Ain -rhs rhs -solu solu -noshift -mat_view
12: with the datafiles in the followig format:
13: Ain (I and J start at 0):
14: ------------------------
15: 3 3 6
16: 0 0 1.0
17: 0 1 2.0
18: 1 0 3.0
19: 1 1 4.0
20: 1 2 5.0
21: 2 2 6.0
23: rhs
24: ---
25: 0 3.0
26: 1 12.0
27: 2 6.0
29: solu
30: ----
31: 0 1.0
32: 0 1.0
33: 0 1.0
34: */
36: #include <petscmat.h>
38: int main(int argc,char **args)
39: {
40: Mat A;
41: Vec b,u,u_tmp;
42: char Ain[PETSC_MAX_PATH_LEN],rhs[PETSC_MAX_PATH_LEN],solu[PETSC_MAX_PATH_LEN];
44: int m,n,nz,dummy; /* these are fscaned so kept as int */
45: PetscInt i,col,row,shift = 1,sizes[3],nsizes;
46: PetscScalar val;
47: PetscReal res_norm;
48: FILE *Afile,*bfile,*ufile;
49: PetscViewer view;
50: PetscBool flg_A,flg_b,flg_u,flg;
51: PetscMPIInt size;
53: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
54: MPI_Comm_size(PETSC_COMM_WORLD,&size);
55: if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!");
57: /* Read in matrix, rhs and exact solution from ascii files */
58: PetscOptionsGetString(NULL,NULL,"-Ain",Ain,PETSC_MAX_PATH_LEN,&flg_A);
59: PetscOptionsHasName(NULL,NULL,"-noshift",&flg);
60: if (flg) shift = 0;
61: if (flg_A) {
62: PetscPrintf(PETSC_COMM_SELF,"\n Read matrix in ascii format ...\n");
63: PetscFOpen(PETSC_COMM_SELF,Ain,"r",&Afile);
64: nsizes = 3;
65: PetscOptionsGetIntArray(NULL,NULL,"-nosizesinfile",sizes,&nsizes,&flg);
66: if (flg) {
67: if (nsizes != 3) SETERRQ(PETSC_COMM_WORLD,1,"Must pass in three m,n,nz as arguments for -nosizesinfile");
68: m = sizes[0];
69: n = sizes[1];
70: nz = sizes[2];
71: } else {
72: if (fscanf(Afile,"%d %d %d\n",&m,&n,&nz) != 3) SETERRQ(PETSC_COMM_SELF,1,"Badly formatted input file\n");
73: }
74: PetscPrintf(PETSC_COMM_SELF,"m: %d, n: %d, nz: %d \n", m,n,nz);
75: if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "Number of rows, cols must be same for this example\n");
76: MatCreate(PETSC_COMM_SELF,&A);
77: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);
78: MatSetFromOptions(A);
79: MatSeqAIJSetPreallocation(A,nz/m,NULL);
80: MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
82: for (i=0; i<nz; i++) {
83: if (fscanf(Afile,"%d %d %le\n",&row,&col,(double*)&val) != 3) SETERRQ(PETSC_COMM_SELF,1,"Badly formatted input file\n");
84: row -= shift; col -= shift; /* set index set starts at 0 */
85: MatSetValues(A,1,&row,1,&col,&val,INSERT_VALUES);
86: }
87: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
88: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
89: fflush(stdout);
90: fclose(Afile);
91: }
93: PetscOptionsGetString(NULL,NULL,"-rhs",rhs,PETSC_MAX_PATH_LEN,&flg_b);
94: if (flg_b) {
95: VecCreate(PETSC_COMM_SELF,&b);
96: VecSetSizes(b,PETSC_DECIDE,n);
97: VecSetFromOptions(b);
98: PetscPrintf(PETSC_COMM_SELF,"\n Read rhs in ascii format ...\n");
99: PetscFOpen(PETSC_COMM_SELF,rhs,"r",&bfile);
100: for (i=0; i<n; i++) {
101: if (fscanf(bfile,"%d %le\n",&dummy,(double*)&val) != 2) SETERRQ(PETSC_COMM_SELF,1,"Badly formatted input file\n");
102: VecSetValues(b,1,&i,&val,INSERT_VALUES);
103: }
104: VecAssemblyBegin(b);
105: VecAssemblyEnd(b);
106: fflush(stdout);
107: fclose(bfile);
108: }
110: PetscOptionsGetString(NULL,NULL,"-solu",solu,PETSC_MAX_PATH_LEN,&flg_u);
111: if (flg_u) {
112: VecCreate(PETSC_COMM_SELF,&u);
113: VecSetSizes(u,PETSC_DECIDE,n);
114: VecSetFromOptions(u);
115: PetscPrintf(PETSC_COMM_SELF,"\n Read exact solution in ascii format ...\n");
116: PetscFOpen(PETSC_COMM_SELF,solu,"r",&ufile);
117: for (i=0; i<n; i++) {
118: if (fscanf(ufile,"%d %le\n",&dummy,(double*)&val) != 2) SETERRQ(PETSC_COMM_SELF,1,"Badly formatted input file\n");
119: VecSetValues(u,1,&i,&val,INSERT_VALUES);
120: }
121: VecAssemblyBegin(u);
122: VecAssemblyEnd(u);
123: fflush(stdout);
124: fclose(ufile);
125: }
127: /* Write matrix, rhs and exact solution in Petsc binary file */
128: PetscPrintf(PETSC_COMM_SELF,"\n Write matrix in binary to 'matrix.dat' ...\n");
129: PetscViewerBinaryOpen(PETSC_COMM_SELF,"matrix.dat",FILE_MODE_WRITE,&view);
130: MatView(A,view);
132: if (flg_b) { /* Write rhs in Petsc binary file */
133: PetscPrintf(PETSC_COMM_SELF,"\n Write rhs in binary to 'matrix.dat' ...\n");
134: VecView(b,view);
135: }
136: if (flg_u) {
137: PetscPrintf(PETSC_COMM_SELF,"\n Write exact solution in binary to 'matrix.dat' ...\n");
138: VecView(u,view);
139: }
141: /* Check accuracy of the data */
142: if (flg_A & flg_b & flg_u) {
143: VecDuplicate(u,&u_tmp);
144: MatMult(A,u,u_tmp);
145: VecAXPY(u_tmp,-1.0,b);
146: VecNorm(u_tmp,NORM_2,&res_norm);
147: PetscPrintf(PETSC_COMM_SELF,"\n Accuracy of the reading data: | b - A*u |_2 : %g \n",res_norm);
148: VecDestroy(&u_tmp);
149: }
151: MatDestroy(&A);
152: if (flg_b) {VecDestroy(&b);}
153: if (flg_u) {VecDestroy(&u);}
154: PetscViewerDestroy(&view);
155: PetscFinalize();
156: return ierr;
157: }