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 */ 179371c9d4SSatish Balay PetscErrorCode PadMatrix(Mat A, Vec v, PetscScalar c, Mat *B) { 18c4762a1bSJed Brown PetscInt n, i, *cnt, *indices, nc; 19c4762a1bSJed Brown const PetscInt *aj; 20c4762a1bSJed Brown const PetscScalar *vv, *aa; 21c4762a1bSJed Brown 22c4762a1bSJed Brown PetscFunctionBegin; 239566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &n, NULL)); 249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(v, &vv)); 259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &indices)); 26c4762a1bSJed Brown for (i = 0; i < n; i++) indices[i] = i; 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* determine number of nonzeros per row in the new matrix */ 299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n + 1, &cnt)); 30c4762a1bSJed Brown for (i = 0; i < n; i++) { 319566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i, &nc, NULL, NULL)); 32c4762a1bSJed Brown cnt[i] = nc + (vv[i] != 0.0); 339566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i, &nc, NULL, NULL)); 34c4762a1bSJed Brown } 35c4762a1bSJed Brown cnt[n] = 1; 36*ad540459SPierre Jolivet for (i = 0; i < n; i++) cnt[n] += (vv[i] != 0.0); 379566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, n + 1, n + 1, 0, cnt, B)); 389566063dSJacob Faibussowitsch PetscCall(MatSetOption(*B, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* copy over the matrix entries from the matrix and then the vector */ 41c4762a1bSJed Brown for (i = 0; i < n; i++) { 429566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i, &nc, &aj, &aa)); 439566063dSJacob Faibussowitsch PetscCall(MatSetValues(*B, 1, &i, nc, aj, aa, INSERT_VALUES)); 449566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i, &nc, &aj, &aa)); 45c4762a1bSJed Brown } 469566063dSJacob Faibussowitsch PetscCall(MatSetValues(*B, 1, &n, n, indices, vv, INSERT_VALUES)); 479566063dSJacob Faibussowitsch PetscCall(MatSetValues(*B, n, indices, 1, &n, vv, INSERT_VALUES)); 489566063dSJacob Faibussowitsch PetscCall(MatSetValues(*B, 1, &n, 1, &n, &c, INSERT_VALUES)); 49c4762a1bSJed Brown 509566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*B, MAT_FINAL_ASSEMBLY)); 519566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*B, MAT_FINAL_ASSEMBLY)); 529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(v, &vv)); 539566063dSJacob Faibussowitsch PetscCall(PetscFree(cnt)); 549566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 55c4762a1bSJed Brown PetscFunctionReturn(0); 56c4762a1bSJed Brown } 57c4762a1bSJed Brown 589371c9d4SSatish Balay int main(int argc, char **args) { 59c4762a1bSJed Brown Mat A, B; 60c4762a1bSJed Brown PetscViewer fd; /* viewer */ 61c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; /* input file name */ 62c4762a1bSJed Brown PetscBool flg; 63c4762a1bSJed Brown Vec v; 64c4762a1bSJed Brown 65327415f7SBarry Smith PetscFunctionBeginUser; 669566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 67c4762a1bSJed Brown /* 68c4762a1bSJed Brown Determine files from which we read the two linear systems 69c4762a1bSJed Brown (matrix and right-hand-side vector). 70c4762a1bSJed Brown */ 719566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-f0", file, sizeof(file), &flg)); 7228b400f6SJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f0 option"); 73c4762a1bSJed Brown 749566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd)); 75c4762a1bSJed Brown 769566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 779566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATSEQAIJ)); 789566063dSJacob Faibussowitsch PetscCall(MatLoad(A, fd)); 799566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &v)); 809566063dSJacob Faibussowitsch PetscCall(VecLoad(v, fd)); 819566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_SELF)); 829566063dSJacob Faibussowitsch PetscCall(PadMatrix(A, v, 3.0, &B)); 839566063dSJacob Faibussowitsch PetscCall(MatView(B, PETSC_VIEWER_STDOUT_SELF)); 849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 879566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 90b122ec5aSJacob Faibussowitsch return 0; 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 93c4762a1bSJed Brown /*TEST 94c4762a1bSJed Brown 95c4762a1bSJed Brown test: 96c4762a1bSJed Brown args: -f0 ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64 97dfd57a17SPierre Jolivet requires: double !complex !defined(PETSC_USE_64BIT_INDICES) 98c4762a1bSJed Brown 99c4762a1bSJed Brown TEST*/ 100