xref: /petsc/src/mat/tests/ex13.c (revision 98d129c30f3ee9fdddc40fdbc5a989b7be64f888)
1 static char help[] = "Tests copying and ordering uniprocessor row-based sparse matrices.\n\n";
2 
3 #include <petscmat.h>
4 
5 int main(int argc, char **args)
6 {
7   Mat         C, A;
8   PetscInt    i, j, m = 5, n = 5, Ii, J;
9   PetscScalar v;
10   IS          perm, iperm;
11 
12   PetscFunctionBeginUser;
13   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
14   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, m * n, m * n, 5, NULL, &C));
15   /* create the matrix for the five point stencil, YET AGAIN*/
16   for (i = 0; i < m; i++) {
17     for (j = 0; j < n; j++) {
18       v  = -1.0;
19       Ii = j + n * i;
20       if (i > 0) {
21         J = Ii - n;
22         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
23       }
24       if (i < m - 1) {
25         J = Ii + n;
26         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
27       }
28       if (j > 0) {
29         J = Ii - 1;
30         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
31       }
32       if (j < n - 1) {
33         J = Ii + 1;
34         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
35       }
36       v = 4.0;
37       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
38     }
39   }
40   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
41   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
42 
43   PetscCall(MatConvert(C, MATSAME, MAT_INITIAL_MATRIX, &A));
44 
45   PetscCall(MatGetOrdering(A, MATORDERINGND, &perm, &iperm));
46   PetscCall(ISView(perm, PETSC_VIEWER_STDOUT_SELF));
47   PetscCall(ISView(iperm, PETSC_VIEWER_STDOUT_SELF));
48   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_SELF));
49 
50   PetscCall(ISDestroy(&perm));
51   PetscCall(ISDestroy(&iperm));
52   PetscCall(MatDestroy(&C));
53   PetscCall(MatDestroy(&A));
54   PetscCall(PetscFinalize());
55   return 0;
56 }
57 
58 /*TEST
59 
60    test:
61       filter: grep -v " MPI process"
62 
63 TEST*/
64