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 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; 54*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 552c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,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 */ 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-Ain",Ain,sizeof(Ain),&flg_A)); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-noshift",&flg)); 60c4762a1bSJed Brown if (flg) shift = 0; 61c4762a1bSJed Brown if (flg_A) { 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n Read matrix in ascii format ...\n")); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFOpen(PETSC_COMM_SELF,Ain,"r",&Afile)); 64c4762a1bSJed Brown nsizes = 3; 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetIntArray(NULL,NULL,"-nosizesinfile",sizes,&nsizes,&flg)); 66c4762a1bSJed Brown if (flg) { 672c71b3e2SJacob Faibussowitsch PetscCheckFalse(nsizes != 3,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 { 722c71b3e2SJacob Faibussowitsch PetscCheckFalse(fscanf(Afile,"%d %d %d\n",&m,&n,&nz) != 3,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Badly formatted input file"); 73c4762a1bSJed Brown } 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"m: %d, n: %d, nz: %d \n", m,n,nz)); 752c71b3e2SJacob Faibussowitsch PetscCheckFalse(m != n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "Number of rows, cols must be same for this example"); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_SELF,&A)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(A,nz/m,NULL)); 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown for (i=0; i<nz; i++) { 832c71b3e2SJacob Faibussowitsch PetscCheckFalse(fscanf(Afile,"%d %d %le\n",&row,&col,(double*)&val) != 3,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Badly formatted input file"); 84c4762a1bSJed Brown row -= shift; col -= shift; /* set index set starts at 0 */ 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&row,1,&col,&val,INSERT_VALUES)); 86c4762a1bSJed Brown } 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 89c4762a1bSJed Brown fclose(Afile); 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-rhs",rhs,sizeof(rhs),&flg_b)); 93c4762a1bSJed Brown if (flg_b) { 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_SELF,&b)); 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(b,PETSC_DECIDE,n)); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(b)); 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n Read rhs in ascii format ...\n")); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFOpen(PETSC_COMM_SELF,rhs,"r",&bfile)); 99c4762a1bSJed Brown for (i=0; i<n; i++) { 1002c71b3e2SJacob Faibussowitsch PetscCheckFalse(fscanf(bfile,"%d %le\n",&dummy,(double*)&val) != 2,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Badly formatted input file"); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(b,1,&i,&val,INSERT_VALUES)); 102c4762a1bSJed Brown } 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(b)); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(b)); 105c4762a1bSJed Brown fclose(bfile); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-solu",solu,sizeof(solu),&flg_u)); 109c4762a1bSJed Brown if (flg_u) { 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_SELF,&u)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(u,PETSC_DECIDE,n)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(u)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n Read exact solution in ascii format ...\n")); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFOpen(PETSC_COMM_SELF,solu,"r",&ufile)); 115c4762a1bSJed Brown for (i=0; i<n; i++) { 1162c71b3e2SJacob Faibussowitsch PetscCheckFalse(fscanf(ufile,"%d %le\n",&dummy,(double*)&val) != 2,PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Badly formatted input file"); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(u,1,&i,&val,INSERT_VALUES)); 118c4762a1bSJed Brown } 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(u)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(u)); 121c4762a1bSJed Brown fclose(ufile); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown 124c4762a1bSJed Brown /* Write matrix, rhs and exact solution in Petsc binary file */ 125*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n Write matrix in binary to 'matrix.dat' ...\n")); 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_SELF,"matrix.dat",FILE_MODE_WRITE,&view)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,view)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown if (flg_b) { /* Write rhs in Petsc binary file */ 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n Write rhs in binary to 'matrix.dat' ...\n")); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(b,view)); 132c4762a1bSJed Brown } 133c4762a1bSJed Brown if (flg_u) { 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n Write exact solution in binary to 'matrix.dat' ...\n")); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(u,view)); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Check accuracy of the data */ 139c4762a1bSJed Brown if (flg_A & flg_b & flg_u) { 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(u,&u_tmp)); 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,u,u_tmp)); 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(u_tmp,-1.0,b)); 143*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(u_tmp,NORM_2,&res_norm)); 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"\n Accuracy of the reading data: | b - A*u |_2 : %g \n",res_norm)); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u_tmp)); 146c4762a1bSJed Brown } 147c4762a1bSJed Brown 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 149*5f80ce2aSJacob Faibussowitsch if (flg_b) CHKERRQ(VecDestroy(&b)); 150*5f80ce2aSJacob Faibussowitsch if (flg_u) CHKERRQ(VecDestroy(&u)); 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 152c4762a1bSJed Brown ierr = PetscFinalize(); 153c4762a1bSJed Brown return ierr; 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156c4762a1bSJed Brown /*TEST 157c4762a1bSJed Brown 158c4762a1bSJed Brown build: 159dfd57a17SPierre 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