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 12a5b23f4aSJose E. Roman with the datafiles in the following 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 int m,n = 0,nz,dummy; /* these are fscaned so kept as int */ 44c4762a1bSJed Brown PetscInt i,col,row,shift = 1,sizes[3],nsizes; 45c4762a1bSJed Brown PetscScalar val; 46c4762a1bSJed Brown PetscReal res_norm; 47c4762a1bSJed Brown FILE *Afile,*bfile,*ufile; 48c4762a1bSJed Brown PetscViewer view; 49c4762a1bSJed Brown PetscBool flg_A,flg_b,flg_u,flg; 50c4762a1bSJed Brown PetscMPIInt size; 51c4762a1bSJed Brown 529566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 54be096a46SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Read in matrix, rhs and exact solution from ascii files */ 579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-Ain",Ain,sizeof(Ain),&flg_A)); 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-noshift",&flg)); 59c4762a1bSJed Brown if (flg) shift = 0; 60c4762a1bSJed Brown if (flg_A) { 619566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n Read matrix in ascii format ...\n")); 629566063dSJacob Faibussowitsch PetscCall(PetscFOpen(PETSC_COMM_SELF,Ain,"r",&Afile)); 63c4762a1bSJed Brown nsizes = 3; 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL,NULL,"-nosizesinfile",sizes,&nsizes,&flg)); 65c4762a1bSJed Brown if (flg) { 66*08401ef6SPierre Jolivet PetscCheck(nsizes == 3,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must pass in three m,n,nz as arguments for -nosizesinfile"); 67c4762a1bSJed Brown m = sizes[0]; 68c4762a1bSJed Brown n = sizes[1]; 69c4762a1bSJed Brown nz = sizes[2]; 70c4762a1bSJed Brown } else { 71*08401ef6SPierre Jolivet PetscCheck(fscanf(Afile,"%d %d %d\n",&m,&n,&nz) == 3,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Badly formatted input file"); 72c4762a1bSJed Brown } 739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"m: %d, n: %d, nz: %d \n", m,n,nz)); 74*08401ef6SPierre Jolivet PetscCheck(m == n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "Number of rows, cols must be same for this example"); 759566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&A)); 769566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n)); 779566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 789566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,nz/m,NULL)); 799566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown for (i=0; i<nz; i++) { 822c71b3e2SJacob Faibussowitsch PetscCheckFalse(fscanf(Afile,"%d %d %le\n",&row,&col,(double*)&val) != 3,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Badly formatted input file"); 83c4762a1bSJed Brown row -= shift; col -= shift; /* set index set starts at 0 */ 849566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&row,1,&col,&val,INSERT_VALUES)); 85c4762a1bSJed Brown } 869566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 879566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 88c4762a1bSJed Brown fclose(Afile); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-rhs",rhs,sizeof(rhs),&flg_b)); 92c4762a1bSJed Brown if (flg_b) { 939566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,&b)); 949566063dSJacob Faibussowitsch PetscCall(VecSetSizes(b,PETSC_DECIDE,n)); 959566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(b)); 969566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n Read rhs in ascii format ...\n")); 979566063dSJacob Faibussowitsch PetscCall(PetscFOpen(PETSC_COMM_SELF,rhs,"r",&bfile)); 98c4762a1bSJed Brown for (i=0; i<n; i++) { 992c71b3e2SJacob Faibussowitsch PetscCheckFalse(fscanf(bfile,"%d %le\n",&dummy,(double*)&val) != 2,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Badly formatted input file"); 1009566063dSJacob Faibussowitsch PetscCall(VecSetValues(b,1,&i,&val,INSERT_VALUES)); 101c4762a1bSJed Brown } 1029566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(b)); 1039566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(b)); 104c4762a1bSJed Brown fclose(bfile); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown 1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-solu",solu,sizeof(solu),&flg_u)); 108c4762a1bSJed Brown if (flg_u) { 1099566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,&u)); 1109566063dSJacob Faibussowitsch PetscCall(VecSetSizes(u,PETSC_DECIDE,n)); 1119566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(u)); 1129566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n Read exact solution in ascii format ...\n")); 1139566063dSJacob Faibussowitsch PetscCall(PetscFOpen(PETSC_COMM_SELF,solu,"r",&ufile)); 114c4762a1bSJed Brown for (i=0; i<n; i++) { 1152c71b3e2SJacob Faibussowitsch PetscCheckFalse(fscanf(ufile,"%d %le\n",&dummy,(double*)&val) != 2,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Badly formatted input file"); 1169566063dSJacob Faibussowitsch PetscCall(VecSetValues(u,1,&i,&val,INSERT_VALUES)); 117c4762a1bSJed Brown } 1189566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(u)); 1199566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(u)); 120c4762a1bSJed Brown fclose(ufile); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* Write matrix, rhs and exact solution in Petsc binary file */ 1249566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n Write matrix in binary to 'matrix.dat' ...\n")); 1259566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF,"matrix.dat",FILE_MODE_WRITE,&view)); 1269566063dSJacob Faibussowitsch PetscCall(MatView(A,view)); 127c4762a1bSJed Brown 128c4762a1bSJed Brown if (flg_b) { /* Write rhs in Petsc binary file */ 1299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n Write rhs in binary to 'matrix.dat' ...\n")); 1309566063dSJacob Faibussowitsch PetscCall(VecView(b,view)); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown if (flg_u) { 1339566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n Write exact solution in binary to 'matrix.dat' ...\n")); 1349566063dSJacob Faibussowitsch PetscCall(VecView(u,view)); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* Check accuracy of the data */ 138c4762a1bSJed Brown if (flg_A & flg_b & flg_u) { 1399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&u_tmp)); 1409566063dSJacob Faibussowitsch PetscCall(MatMult(A,u,u_tmp)); 1419566063dSJacob Faibussowitsch PetscCall(VecAXPY(u_tmp,-1.0,b)); 1429566063dSJacob Faibussowitsch PetscCall(VecNorm(u_tmp,NORM_2,&res_norm)); 1439566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n Accuracy of the reading data: | b - A*u |_2 : %g \n",res_norm)); 1449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u_tmp)); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 1479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1489566063dSJacob Faibussowitsch if (flg_b) PetscCall(VecDestroy(&b)); 1499566063dSJacob Faibussowitsch if (flg_u) PetscCall(VecDestroy(&u)); 1509566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 1519566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 152b122ec5aSJacob Faibussowitsch return 0; 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown /*TEST 156c4762a1bSJed Brown 157c4762a1bSJed Brown build: 158dfd57a17SPierre Jolivet requires: !defined(PETSC_USE_64BIT_INDICES) double !complex 159c4762a1bSJed Brown 160c4762a1bSJed Brown test: 1615a0df56aSBarry Smith requires: datafilespath 162c4762a1bSJed Brown args: -Ain ${DATAFILESPATH}/matrices/indefinite/afiro_A.dat 163c4762a1bSJed Brown 164c4762a1bSJed Brown TEST*/ 165