xref: /petsc/src/mat/tests/ex14.c (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
1 
2 static char help[] = "Tests sequential and parallel MatGetRow() and MatRestoreRow().\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc, char **args)
7 {
8   Mat                C;
9   PetscInt           i, j, m = 5, n = 5, Ii, J, nz, rstart, rend;
10   PetscMPIInt        rank;
11   const PetscInt    *idx;
12   PetscScalar        v;
13   const PetscScalar *values;
14 
15   PetscFunctionBeginUser;
16   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
17   /* Create the matrix for the five point stencil, YET AGAIN */
18   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
19   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
20   PetscCall(MatSetFromOptions(C));
21   PetscCall(MatSetUp(C));
22   for (i = 0; i < m; i++) {
23     for (j = 0; j < n; j++) {
24       v  = -1.0;
25       Ii = j + n * i;
26       if (i > 0) {
27         J = Ii - n;
28         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
29       }
30       if (i < m - 1) {
31         J = Ii + n;
32         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
33       }
34       if (j > 0) {
35         J = Ii - 1;
36         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
37       }
38       if (j < n - 1) {
39         J = Ii + 1;
40         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
41       }
42       v = 4.0;
43       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
44     }
45   }
46   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
47   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
48   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
49 
50   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
51   PetscCall(MatGetOwnershipRange(C, &rstart, &rend));
52   for (i = rstart; i < rend; i++) {
53     PetscCall(MatGetRow(C, i, &nz, &idx, &values));
54     if (rank == 0) {
55       for (j = 0; j < nz; j++) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT " %g ", idx[j], (double)PetscRealPart(values[j])));
56       PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
57     }
58     PetscCall(MatRestoreRow(C, i, &nz, &idx, &values));
59   }
60 
61   PetscCall(MatDestroy(&C));
62   PetscCall(PetscFinalize());
63   return 0;
64 }
65 
66 /*TEST
67 
68    test:
69 
70 TEST*/
71