1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Reads a PETSc matrix and vector from a file; expands the matrix with the vector\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* 5c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices. 6c4762a1bSJed Brown automatically includes: 7c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 8c4762a1bSJed Brown petscmat.h - matrices 9c4762a1bSJed Brown petscis.h - index sets petscviewer.h - viewers 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown #include <petscmat.h> 12c4762a1bSJed Brown 13c4762a1bSJed Brown /* 14c4762a1bSJed Brown 15c4762a1bSJed Brown Adds a new column and row to the vector (the last) containing the vector 16c4762a1bSJed Brown */ 17c4762a1bSJed Brown PetscErrorCode PadMatrix(Mat A,Vec v,PetscScalar c,Mat *B) 18c4762a1bSJed Brown { 19c4762a1bSJed Brown PetscInt n,i,*cnt,*indices,nc; 20c4762a1bSJed Brown const PetscInt *aj; 21c4762a1bSJed Brown const PetscScalar *vv,*aa; 22c4762a1bSJed Brown 23c4762a1bSJed Brown PetscFunctionBegin; 249566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&n,NULL)); 259566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(v,&vv)); 269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&indices)); 27c4762a1bSJed Brown for (i=0; i<n; i++) indices[i] = i; 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* determine number of nonzeros per row in the new matrix */ 309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n+1,&cnt)); 31c4762a1bSJed Brown for (i=0; i<n; i++) { 329566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i,&nc,NULL,NULL)); 33c4762a1bSJed Brown cnt[i] = nc + (vv[i] != 0.0); 349566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i,&nc,NULL,NULL)); 35c4762a1bSJed Brown } 36c4762a1bSJed Brown cnt[n] = 1; 37c4762a1bSJed Brown for (i=0; i<n; i++) { 38c4762a1bSJed Brown cnt[n] += (vv[i] != 0.0); 39c4762a1bSJed Brown } 409566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,n+1,n+1,0,cnt,B)); 419566063dSJacob Faibussowitsch PetscCall(MatSetOption(*B,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* copy over the matrix entries from the matrix and then the vector */ 44c4762a1bSJed Brown for (i=0; i<n; i++) { 459566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,i,&nc,&aj,&aa)); 469566063dSJacob Faibussowitsch PetscCall(MatSetValues(*B,1,&i,nc,aj,aa,INSERT_VALUES)); 479566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,i,&nc,&aj,&aa)); 48c4762a1bSJed Brown } 499566063dSJacob Faibussowitsch PetscCall(MatSetValues(*B,1,&n,n,indices,vv,INSERT_VALUES)); 509566063dSJacob Faibussowitsch PetscCall(MatSetValues(*B,n,indices,1,&n,vv,INSERT_VALUES)); 519566063dSJacob Faibussowitsch PetscCall(MatSetValues(*B,1,&n,1,&n,&c,INSERT_VALUES)); 52c4762a1bSJed Brown 539566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY)); 549566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY)); 559566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(v,&vv)); 569566063dSJacob Faibussowitsch PetscCall(PetscFree(cnt)); 579566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 58c4762a1bSJed Brown PetscFunctionReturn(0); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61c4762a1bSJed Brown int main(int argc,char **args) 62c4762a1bSJed Brown { 63c4762a1bSJed Brown Mat A,B; 64c4762a1bSJed Brown PetscViewer fd; /* viewer */ 65c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; /* input file name */ 66c4762a1bSJed Brown PetscBool flg; 67c4762a1bSJed Brown Vec v; 68c4762a1bSJed Brown 69*327415f7SBarry Smith PetscFunctionBeginUser; 709566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 71c4762a1bSJed Brown /* 72c4762a1bSJed Brown Determine files from which we read the two linear systems 73c4762a1bSJed Brown (matrix and right-hand-side vector). 74c4762a1bSJed Brown */ 759566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-f0",file,sizeof(file),&flg)); 7628b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f0 option"); 77c4762a1bSJed Brown 789566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 79c4762a1bSJed Brown 809566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 819566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSEQAIJ)); 829566063dSJacob Faibussowitsch PetscCall(MatLoad(A,fd)); 839566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&v)); 849566063dSJacob Faibussowitsch PetscCall(VecLoad(v,fd)); 859566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF)); 869566063dSJacob Faibussowitsch PetscCall(PadMatrix(A,v,3.0,&B)); 879566063dSJacob Faibussowitsch PetscCall(MatView(B,PETSC_VIEWER_STDOUT_SELF)); 889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 919566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 92c4762a1bSJed Brown 939566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 94b122ec5aSJacob Faibussowitsch return 0; 95c4762a1bSJed Brown } 96c4762a1bSJed Brown 97c4762a1bSJed Brown /*TEST 98c4762a1bSJed Brown 99c4762a1bSJed Brown test: 100c4762a1bSJed Brown args: -f0 ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64 101dfd57a17SPierre Jolivet requires: double !complex !defined(PETSC_USE_64BIT_INDICES) 102c4762a1bSJed Brown 103c4762a1bSJed Brown TEST*/ 104