xref: /petsc/src/mat/tests/ex14.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests sequential and parallel MatGetRow() and MatRestoreRow().\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6*9371c9d4SSatish Balay int main(int argc, char **args) {
7c4762a1bSJed Brown   Mat                C;
8c4762a1bSJed Brown   PetscInt           i, j, m = 5, n = 5, Ii, J, nz, rstart, rend;
9c4762a1bSJed Brown   PetscMPIInt        rank;
10c4762a1bSJed Brown   const PetscInt    *idx;
11c4762a1bSJed Brown   PetscScalar        v;
12c4762a1bSJed Brown   const PetscScalar *values;
13c4762a1bSJed Brown 
14327415f7SBarry Smith   PetscFunctionBeginUser;
159566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
16c4762a1bSJed Brown   /* Create the matrix for the five point stencil, YET AGAIN */
179566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
189566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
199566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
209566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
21c4762a1bSJed Brown   for (i = 0; i < m; i++) {
22c4762a1bSJed Brown     for (j = 0; j < n; j++) {
23*9371c9d4SSatish Balay       v  = -1.0;
24*9371c9d4SSatish Balay       Ii = j + n * i;
25*9371c9d4SSatish Balay       if (i > 0) {
26*9371c9d4SSatish Balay         J = Ii - n;
27*9371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
28*9371c9d4SSatish Balay       }
29*9371c9d4SSatish Balay       if (i < m - 1) {
30*9371c9d4SSatish Balay         J = Ii + n;
31*9371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
32*9371c9d4SSatish Balay       }
33*9371c9d4SSatish Balay       if (j > 0) {
34*9371c9d4SSatish Balay         J = Ii - 1;
35*9371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
36*9371c9d4SSatish Balay       }
37*9371c9d4SSatish Balay       if (j < n - 1) {
38*9371c9d4SSatish Balay         J = Ii + 1;
39*9371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
40*9371c9d4SSatish Balay       }
41*9371c9d4SSatish Balay       v = 4.0;
42*9371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
43c4762a1bSJed Brown     }
44c4762a1bSJed Brown   }
459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
479566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
48c4762a1bSJed Brown 
499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
509566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(C, &rstart, &rend));
51c4762a1bSJed Brown   for (i = rstart; i < rend; i++) {
529566063dSJacob Faibussowitsch     PetscCall(MatGetRow(C, i, &nz, &idx, &values));
53dd400576SPatrick Sanan     if (rank == 0) {
54*9371c9d4SSatish Balay       for (j = 0; j < nz; j++) { PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT " %g ", idx[j], (double)PetscRealPart(values[j]))); }
559566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
56c4762a1bSJed Brown     }
579566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(C, i, &nz, &idx, &values));
58c4762a1bSJed Brown   }
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
619566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
62b122ec5aSJacob Faibussowitsch   return 0;
63c4762a1bSJed Brown }
64c4762a1bSJed Brown 
65c4762a1bSJed Brown /*TEST
66c4762a1bSJed Brown 
67c4762a1bSJed Brown    test:
68c4762a1bSJed Brown 
69c4762a1bSJed Brown TEST*/
70