xref: /petsc/src/mat/tests/ex10.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests repeated use of assembly for matrices.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6*9371c9d4SSatish Balay int main(int argc, char **args) {
7c4762a1bSJed Brown   Mat         C;
8c4762a1bSJed Brown   PetscInt    i, j, m = 5, n = 2, Ii, J;
9c4762a1bSJed Brown   PetscMPIInt rank, size;
10c4762a1bSJed Brown   PetscScalar v;
11c4762a1bSJed Brown 
12327415f7SBarry Smith   PetscFunctionBeginUser;
139566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
16c4762a1bSJed Brown   n = 2 * size;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN*/
199566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
209566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
219566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
229566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
23c4762a1bSJed Brown   for (i = 0; i < m; i++) {
24c4762a1bSJed Brown     for (j = 2 * rank; j < 2 * rank + 2; j++) {
25*9371c9d4SSatish Balay       v  = -1.0;
26*9371c9d4SSatish Balay       Ii = j + n * i;
27*9371c9d4SSatish Balay       if (i > 0) {
28*9371c9d4SSatish Balay         J = Ii - n;
29*9371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
30*9371c9d4SSatish Balay       }
31*9371c9d4SSatish Balay       if (i < m - 1) {
32*9371c9d4SSatish Balay         J = Ii + n;
33*9371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
34*9371c9d4SSatish Balay       }
35*9371c9d4SSatish Balay       if (j > 0) {
36*9371c9d4SSatish Balay         J = Ii - 1;
37*9371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
38*9371c9d4SSatish Balay       }
39*9371c9d4SSatish Balay       if (j < n - 1) {
40*9371c9d4SSatish Balay         J = Ii + 1;
41*9371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
42*9371c9d4SSatish Balay       }
43*9371c9d4SSatish Balay       v = 4.0;
44*9371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
45c4762a1bSJed Brown     }
46c4762a1bSJed Brown   }
479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
49c4762a1bSJed Brown   for (i = 0; i < m; i++) {
50c4762a1bSJed Brown     for (j = 2 * rank; j < 2 * rank + 2; j++) {
51*9371c9d4SSatish Balay       v  = 1.0;
52*9371c9d4SSatish Balay       Ii = j + n * i;
53*9371c9d4SSatish Balay       if (i > 0) {
54*9371c9d4SSatish Balay         J = Ii - n;
55*9371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
56*9371c9d4SSatish Balay       }
57*9371c9d4SSatish Balay       if (i < m - 1) {
58*9371c9d4SSatish Balay         J = Ii + n;
59*9371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
60*9371c9d4SSatish Balay       }
61*9371c9d4SSatish Balay       if (j > 0) {
62*9371c9d4SSatish Balay         J = Ii - 1;
63*9371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
64*9371c9d4SSatish Balay       }
65*9371c9d4SSatish Balay       if (j < n - 1) {
66*9371c9d4SSatish Balay         J = Ii + 1;
67*9371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
68*9371c9d4SSatish Balay       }
69*9371c9d4SSatish Balay       v = -4.0;
70*9371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
71c4762a1bSJed Brown     }
72c4762a1bSJed Brown   }
739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
75c4762a1bSJed Brown 
769566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
799566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
80b122ec5aSJacob Faibussowitsch   return 0;
81c4762a1bSJed Brown }
82c4762a1bSJed Brown 
83c4762a1bSJed Brown /*TEST
84c4762a1bSJed Brown 
85c4762a1bSJed Brown    test:
86c4762a1bSJed Brown 
87c4762a1bSJed Brown TEST*/
88