Actual source code: ex50.c
petsc-3.4.5 2014-06-29
2: static char help[] = "Reads in a matrix and vector in ASCII format. Writes\n\
3: them using the PETSc sparse format. Input parameters are:\n\
4: -fin <filename> : input file\n\
5: -fout <filename> : output file\n\n";
7: #include <petscmat.h>
11: int main(int argc,char **args)
12: {
13: Mat A;
14: Vec b;
15: char filein[PETSC_MAX_PATH_LEN],finname[PETSC_MAX_PATH_LEN],fileout[PETSC_MAX_PATH_LEN];
16: PetscInt n,col,row,rowin;
18: PetscBool flg;
19: PetscScalar val,*array;
20: FILE *file;
21: PetscViewer view;
23: PetscInitialize(&argc,&args,(char*)0,help);
25: /* Read in matrix and RHS */
26: PetscOptionsGetString(NULL,"-fin",filein,256,&flg);
27: if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Must indicate file for reading");
28: PetscOptionsGetString(NULL,"-fout",fileout,256,&flg);
29: if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Must indicate file for writing");
31: PetscFixFilename(filein,finname);
32: if (!(file = fopen(finname,"r"))) SETERRQ(PETSC_COMM_SELF,1,"cannot open input file\n");
33: fscanf(file,"%d\n",&n);
35: MatCreate(PETSC_COMM_WORLD,&A);
36: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
37: MatSetFromOptions(A);
38: VecCreate(PETSC_COMM_WORLD,&b);
39: VecSetSizes(b,PETSC_DECIDE,n);
40: VecSetFromOptions(b);
42: for (row=0; row<n; row++) {
43: fscanf(file,"row %d:",&rowin);
44: if (rowin != row) SETERRQ(PETSC_COMM_SELF,1,"Bad file");
45: while (fscanf(file," %d %le",&col,(double*)&val)) {
46: MatSetValues(A,1,&row,1,&col,&val,INSERT_VALUES);
47: }
48: }
49: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
50: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
51: VecGetArray(b,&array);
52: for (row=0; row<n; row++) {
53: fscanf(file," ii= %d %le",&col,(double*)(array+row));
54: }
55: VecRestoreArray(b,&array);
57: fclose(file);
59: PetscPrintf(PETSC_COMM_SELF,"Reading matrix complete.\n");
60: PetscViewerBinaryOpen(PETSC_COMM_WORLD,fileout,FILE_MODE_WRITE,&view);
61: MatView(A,view);
62: VecView(b,view);
63: PetscViewerDestroy(&view);
65: VecDestroy(&b);
66: MatDestroy(&A);
68: PetscFinalize();
69: return 0;
70: }