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