xref: /petsc/src/mat/tests/ex8.c (revision 856bee69f0e0908e75ff837867b1777dfb1ced96)
1 
2 static char help[] = "Tests automatic allocation of matrix storage space.\n\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc, char **args)
7 {
8   Mat         C;
9   PetscInt    i, j, m = 3, n = 3, Ii, J;
10   PetscScalar v;
11   MatInfo     info;
12 
13   PetscFunctionBeginUser;
14   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
15   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
16   PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
17 
18   /* create the matrix for the five point stencil, YET AGAIN */
19   PetscCall(MatCreate(PETSC_COMM_SELF, &C));
20   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
21   PetscCall(MatSetFromOptions(C));
22   PetscCall(MatSetUp(C));
23   for (i = 0; i < m; i++) {
24     for (j = 0; j < n; j++) {
25       v  = -1.0;
26       Ii = j + n * i;
27       if (i > 0) {
28         J = Ii - n;
29         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
30       }
31       if (i < m - 1) {
32         J = Ii + n;
33         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
34       }
35       if (j > 0) {
36         J = Ii - 1;
37         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
38       }
39       if (j < n - 1) {
40         J = Ii + 1;
41         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
42       }
43       v = 4.0;
44       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
45     }
46   }
47   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
48   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
49   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_SELF));
50 
51   PetscCall(MatGetInfo(C, MAT_LOCAL, &info));
52   PetscCall(PetscPrintf(PETSC_COMM_SELF, "matrix nonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n", (PetscInt)info.nz_used, (PetscInt)info.nz_allocated));
53 
54   PetscCall(MatDestroy(&C));
55   PetscCall(PetscFinalize());
56   return 0;
57 }
58 
59 /*TEST
60 
61    test:
62 
63 TEST*/
64