xref: /petsc/src/mat/tests/ex10.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 
2 static char help[] = "Tests repeated use of assembly for matrices.\n\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc, char **args) {
7   Mat         C;
8   PetscInt    i, j, m = 5, n = 2, Ii, J;
9   PetscMPIInt rank, size;
10   PetscScalar v;
11 
12   PetscFunctionBeginUser;
13   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
14   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
15   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
16   n = 2 * size;
17 
18   /* create the matrix for the five point stencil, YET AGAIN*/
19   PetscCall(MatCreate(PETSC_COMM_WORLD, &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 = 2 * rank; j < 2 * rank + 2; 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   for (i = 0; i < m; i++) {
50     for (j = 2 * rank; j < 2 * rank + 2; j++) {
51       v  = 1.0;
52       Ii = j + n * i;
53       if (i > 0) {
54         J = Ii - n;
55         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
56       }
57       if (i < m - 1) {
58         J = Ii + n;
59         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
60       }
61       if (j > 0) {
62         J = Ii - 1;
63         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
64       }
65       if (j < n - 1) {
66         J = Ii + 1;
67         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
68       }
69       v = -4.0;
70       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
71     }
72   }
73   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
74   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
75 
76   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
77 
78   PetscCall(MatDestroy(&C));
79   PetscCall(PetscFinalize());
80   return 0;
81 }
82 
83 /*TEST
84 
85    test:
86 
87 TEST*/
88