xref: /petsc/src/mat/tests/ex10.c (revision 856bee69f0e0908e75ff837867b1777dfb1ced96)
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, B;
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   PetscCall(MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, &B)); /* test that SeqAIJ non-preallocated matrices can be duplicated */
25   PetscCall(MatDestroy(&C));
26   C = B;
27   for (i = 0; i < m; i++) {
28     for (j = 2 * rank; j < 2 * rank + 2; j++) {
29       v  = -1.0;
30       Ii = j + n * i;
31       if (i > 0) {
32         J = Ii - n;
33         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
34       }
35       if (i < m - 1) {
36         J = Ii + n;
37         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
38       }
39       if (j > 0) {
40         J = Ii - 1;
41         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
42       }
43       if (j < n - 1) {
44         J = Ii + 1;
45         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
46       }
47       v = 4.0;
48       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
49     }
50   }
51   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
52   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
53   for (i = 0; i < m; i++) {
54     for (j = 2 * rank; j < 2 * rank + 2; j++) {
55       v  = 1.0;
56       Ii = j + n * i;
57       if (i > 0) {
58         J = Ii - n;
59         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
60       }
61       if (i < m - 1) {
62         J = Ii + n;
63         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
64       }
65       if (j > 0) {
66         J = Ii - 1;
67         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
68       }
69       if (j < n - 1) {
70         J = Ii + 1;
71         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
72       }
73       v = -4.0;
74       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
75     }
76   }
77   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
78   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
79 
80   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
81 
82   PetscCall(MatDestroy(&C));
83   PetscCall(PetscFinalize());
84   return 0;
85 }
86 
87 /*TEST
88 
89    test:
90 
91 TEST*/
92