xref: /petsc/src/mat/tests/ex10.c (revision 73fdd05bb67e49f40fd8fd311695ff6fdf0b9b8a)
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 {
8   Mat         C;
9   PetscInt    i, j, m = 5, n = 2, Ii, J;
10   PetscMPIInt rank, size;
11   PetscScalar v;
12 
13   PetscFunctionBeginUser;
14   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
15   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
16   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
17   n = 2 * size;
18 
19   /* create the matrix for the five point stencil, YET AGAIN*/
20   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
21   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
22   PetscCall(MatSetFromOptions(C));
23   PetscCall(MatSetUp(C));
24   for (i = 0; i < m; i++) {
25     for (j = 2 * rank; j < 2 * rank + 2; j++) {
26       v  = -1.0;
27       Ii = j + n * i;
28       if (i > 0) {
29         J = Ii - n;
30         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
31       }
32       if (i < m - 1) {
33         J = Ii + n;
34         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
35       }
36       if (j > 0) {
37         J = Ii - 1;
38         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
39       }
40       if (j < n - 1) {
41         J = Ii + 1;
42         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
43       }
44       v = 4.0;
45       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
46     }
47   }
48   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
49   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
50   for (i = 0; i < m; i++) {
51     for (j = 2 * rank; j < 2 * rank + 2; j++) {
52       v  = 1.0;
53       Ii = j + n * i;
54       if (i > 0) {
55         J = Ii - n;
56         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
57       }
58       if (i < m - 1) {
59         J = Ii + n;
60         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
61       }
62       if (j > 0) {
63         J = Ii - 1;
64         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
65       }
66       if (j < n - 1) {
67         J = Ii + 1;
68         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
69       }
70       v = -4.0;
71       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
72     }
73   }
74   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
75   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
76 
77   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
78 
79   PetscCall(MatDestroy(&C));
80   PetscCall(PetscFinalize());
81   return 0;
82 }
83 
84 /*TEST
85 
86    test:
87 
88 TEST*/
89