xref: /petsc/src/mat/tests/ex14.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7*d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat                C;
9c4762a1bSJed Brown   PetscInt           i, j, m = 5, n = 5, Ii, J, nz, rstart, rend;
10c4762a1bSJed Brown   PetscMPIInt        rank;
11c4762a1bSJed Brown   const PetscInt    *idx;
12c4762a1bSJed Brown   PetscScalar        v;
13c4762a1bSJed Brown   const PetscScalar *values;
14c4762a1bSJed Brown 
15327415f7SBarry Smith   PetscFunctionBeginUser;
169566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
17c4762a1bSJed Brown   /* Create the matrix for the five point stencil, YET AGAIN */
189566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
199566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
209566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
219566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
22c4762a1bSJed Brown   for (i = 0; i < m; i++) {
23c4762a1bSJed Brown     for (j = 0; j < n; j++) {
249371c9d4SSatish Balay       v  = -1.0;
259371c9d4SSatish Balay       Ii = j + n * i;
269371c9d4SSatish Balay       if (i > 0) {
279371c9d4SSatish Balay         J = Ii - n;
289371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
299371c9d4SSatish Balay       }
309371c9d4SSatish Balay       if (i < m - 1) {
319371c9d4SSatish Balay         J = Ii + n;
329371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
339371c9d4SSatish Balay       }
349371c9d4SSatish Balay       if (j > 0) {
359371c9d4SSatish Balay         J = Ii - 1;
369371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
379371c9d4SSatish Balay       }
389371c9d4SSatish Balay       if (j < n - 1) {
399371c9d4SSatish Balay         J = Ii + 1;
409371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
419371c9d4SSatish Balay       }
429371c9d4SSatish Balay       v = 4.0;
439371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
44c4762a1bSJed Brown     }
45c4762a1bSJed Brown   }
469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
489566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
519566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(C, &rstart, &rend));
52c4762a1bSJed Brown   for (i = rstart; i < rend; i++) {
539566063dSJacob Faibussowitsch     PetscCall(MatGetRow(C, i, &nz, &idx, &values));
54dd400576SPatrick Sanan     if (rank == 0) {
5548a46eb9SPierre Jolivet       for (j = 0; j < nz; j++) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT " %g ", idx[j], (double)PetscRealPart(values[j])));
569566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
57c4762a1bSJed Brown     }
589566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(C, i, &nz, &idx, &values));
59c4762a1bSJed Brown   }
60c4762a1bSJed Brown 
619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
629566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
63b122ec5aSJacob Faibussowitsch   return 0;
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown /*TEST
67c4762a1bSJed Brown 
68c4762a1bSJed Brown    test:
69c4762a1bSJed Brown 
70c4762a1bSJed Brown TEST*/
71