xref: /petsc/src/mat/tests/ex10.c (revision 31fe6a7d0c059274179307fb538b31c7bd22b9f4)
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 
6d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
7d71ae5a4SJacob Faibussowitsch {
8*31fe6a7dSBarry Smith   Mat         C, B;
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));
24*31fe6a7dSBarry Smith   PetscCall(MatDuplicate(C, MAT_DO_NOT_COPY_VALUES, &B)); /* test that SeqAIJ non-preallocated matrices can be duplicated */
25*31fe6a7dSBarry Smith   PetscCall(MatDestroy(&C));
26*31fe6a7dSBarry Smith   C = B;
27c4762a1bSJed Brown   for (i = 0; i < m; i++) {
28c4762a1bSJed Brown     for (j = 2 * rank; j < 2 * rank + 2; j++) {
299371c9d4SSatish Balay       v  = -1.0;
309371c9d4SSatish Balay       Ii = j + n * i;
319371c9d4SSatish Balay       if (i > 0) {
329371c9d4SSatish Balay         J = Ii - n;
339371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
349371c9d4SSatish Balay       }
359371c9d4SSatish Balay       if (i < m - 1) {
369371c9d4SSatish Balay         J = Ii + n;
379371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
389371c9d4SSatish Balay       }
399371c9d4SSatish Balay       if (j > 0) {
409371c9d4SSatish Balay         J = Ii - 1;
419371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
429371c9d4SSatish Balay       }
439371c9d4SSatish Balay       if (j < n - 1) {
449371c9d4SSatish Balay         J = Ii + 1;
459371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
469371c9d4SSatish Balay       }
479371c9d4SSatish Balay       v = 4.0;
489371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
49c4762a1bSJed Brown     }
50c4762a1bSJed Brown   }
519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
53c4762a1bSJed Brown   for (i = 0; i < m; i++) {
54c4762a1bSJed Brown     for (j = 2 * rank; j < 2 * rank + 2; j++) {
559371c9d4SSatish Balay       v  = 1.0;
569371c9d4SSatish Balay       Ii = j + n * i;
579371c9d4SSatish Balay       if (i > 0) {
589371c9d4SSatish Balay         J = Ii - n;
599371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
609371c9d4SSatish Balay       }
619371c9d4SSatish Balay       if (i < m - 1) {
629371c9d4SSatish Balay         J = Ii + n;
639371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
649371c9d4SSatish Balay       }
659371c9d4SSatish Balay       if (j > 0) {
669371c9d4SSatish Balay         J = Ii - 1;
679371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
689371c9d4SSatish Balay       }
699371c9d4SSatish Balay       if (j < n - 1) {
709371c9d4SSatish Balay         J = Ii + 1;
719371c9d4SSatish Balay         PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
729371c9d4SSatish Balay       }
739371c9d4SSatish Balay       v = -4.0;
749371c9d4SSatish Balay       PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
75c4762a1bSJed Brown     }
76c4762a1bSJed Brown   }
779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
789566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
79c4762a1bSJed Brown 
809566063dSJacob Faibussowitsch   PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
81c4762a1bSJed Brown 
829566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
839566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
84b122ec5aSJacob Faibussowitsch   return 0;
85c4762a1bSJed Brown }
86c4762a1bSJed Brown 
87c4762a1bSJed Brown /*TEST
88c4762a1bSJed Brown 
89c4762a1bSJed Brown    test:
90c4762a1bSJed Brown 
91c4762a1bSJed Brown TEST*/
92