xref: /petsc/src/mat/tests/ex10.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7*d71ae5a4SJacob Faibussowitsch {
8c4762a1bSJed Brown   Mat         C;
9c4762a1bSJed Brown   PetscInt    i, j, m = 5, n = 2, Ii, J;
10c4762a1bSJed Brown   PetscMPIInt rank, size;
11c4762a1bSJed Brown   PetscScalar v;
12c4762a1bSJed Brown 
13327415f7SBarry Smith   PetscFunctionBeginUser;
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
17c4762a1bSJed Brown   n = 2 * size;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   /* create the matrix for the five point stencil, YET AGAIN*/
209566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
219566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
229566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
239566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
24c4762a1bSJed Brown   for (i = 0; i < m; i++) {
25c4762a1bSJed Brown     for (j = 2 * rank; j < 2 * rank + 2; j++) {
269371c9d4SSatish Balay       v  = -1.0;
279371c9d4SSatish Balay       Ii = j + n * i;
289371c9d4SSatish Balay       if (i > 0) {
299371c9d4SSatish Balay         J = Ii - n;
309371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
319371c9d4SSatish Balay       }
329371c9d4SSatish Balay       if (i < m - 1) {
339371c9d4SSatish Balay         J = Ii + n;
349371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
359371c9d4SSatish Balay       }
369371c9d4SSatish Balay       if (j > 0) {
379371c9d4SSatish Balay         J = Ii - 1;
389371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
399371c9d4SSatish Balay       }
409371c9d4SSatish Balay       if (j < n - 1) {
419371c9d4SSatish Balay         J = Ii + 1;
429371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
439371c9d4SSatish Balay       }
449371c9d4SSatish Balay       v = 4.0;
459371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
46c4762a1bSJed Brown     }
47c4762a1bSJed Brown   }
489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
50c4762a1bSJed Brown   for (i = 0; i < m; i++) {
51c4762a1bSJed Brown     for (j = 2 * rank; j < 2 * rank + 2; j++) {
529371c9d4SSatish Balay       v  = 1.0;
539371c9d4SSatish Balay       Ii = j + n * i;
549371c9d4SSatish Balay       if (i > 0) {
559371c9d4SSatish Balay         J = Ii - n;
569371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
579371c9d4SSatish Balay       }
589371c9d4SSatish Balay       if (i < m - 1) {
599371c9d4SSatish Balay         J = Ii + n;
609371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
619371c9d4SSatish Balay       }
629371c9d4SSatish Balay       if (j > 0) {
639371c9d4SSatish Balay         J = Ii - 1;
649371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
659371c9d4SSatish Balay       }
669371c9d4SSatish Balay       if (j < n - 1) {
679371c9d4SSatish Balay         J = Ii + 1;
689371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
699371c9d4SSatish Balay       }
709371c9d4SSatish Balay       v = -4.0;
719371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
72c4762a1bSJed Brown     }
73c4762a1bSJed Brown   }
749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
809566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
81b122ec5aSJacob Faibussowitsch   return 0;
82c4762a1bSJed Brown }
83c4762a1bSJed Brown 
84c4762a1bSJed Brown /*TEST
85c4762a1bSJed Brown 
86c4762a1bSJed Brown    test:
87c4762a1bSJed Brown 
88c4762a1bSJed Brown TEST*/
89