xref: /petsc/src/mat/tests/ex86.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "Testing MatCreateMPIMatConcatenateSeqMat().\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
main(int argc,char ** argv)4d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
5d71ae5a4SJacob Faibussowitsch {
6c4762a1bSJed Brown   Mat         seqmat, mpimat;
7c4762a1bSJed Brown   PetscMPIInt rank;
8c4762a1bSJed Brown   PetscScalar value[3], *vals;
9c4762a1bSJed Brown   PetscInt    i, col[3], n = 5, bs = 1;
10c4762a1bSJed Brown 
11327415f7SBarry Smith   PetscFunctionBeginUser;
12*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   /* Create seqaij matrices of size (n+rank) by n */
179566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &seqmat));
189566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(seqmat, (n + rank) * bs, PETSC_DECIDE, PETSC_DECIDE, n * bs));
199566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(seqmat));
209566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(seqmat, 3 * bs, NULL));
219566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(seqmat, bs, 3, NULL));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   if (bs == 1) {
249371c9d4SSatish Balay     value[0] = -1.0;
259371c9d4SSatish Balay     value[1] = 2.0;
269371c9d4SSatish Balay     value[2] = -1.0;
27c4762a1bSJed Brown     for (i = 1; i < n - 1; i++) {
289371c9d4SSatish Balay       col[0] = i - 1;
299371c9d4SSatish Balay       col[1] = i;
309371c9d4SSatish Balay       col[2] = i + 1;
319566063dSJacob Faibussowitsch       PetscCall(MatSetValues(seqmat, 1, &i, 3, col, value, INSERT_VALUES));
32c4762a1bSJed Brown     }
339371c9d4SSatish Balay     i      = n - 1;
349371c9d4SSatish Balay     col[0] = n - 2;
359371c9d4SSatish Balay     col[1] = n - 1;
369566063dSJacob Faibussowitsch     PetscCall(MatSetValues(seqmat, 1, &i, 2, col, value, INSERT_VALUES));
37c4762a1bSJed Brown 
389371c9d4SSatish Balay     i        = 0;
399371c9d4SSatish Balay     col[0]   = 0;
409371c9d4SSatish Balay     col[1]   = 1;
419371c9d4SSatish Balay     value[0] = 2.0;
429371c9d4SSatish Balay     value[1] = -1.0;
439566063dSJacob Faibussowitsch     PetscCall(MatSetValues(seqmat, 1, &i, 2, col, value, INSERT_VALUES));
44c4762a1bSJed Brown   } else {
45c4762a1bSJed Brown     PetscInt *rows, *cols, j;
469566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(bs * bs, &vals, bs, &rows, bs, &cols));
47c4762a1bSJed Brown     /* diagonal blocks */
48c4762a1bSJed Brown     for (i = 0; i < bs * bs; i++) vals[i] = 2.0;
49c4762a1bSJed Brown     for (i = 0; i < n * bs; i += bs) {
509371c9d4SSatish Balay       for (j = 0; j < bs; j++) {
519371c9d4SSatish Balay         rows[j] = i + j;
529371c9d4SSatish Balay         cols[j] = i + j;
539371c9d4SSatish Balay       }
549566063dSJacob Faibussowitsch       PetscCall(MatSetValues(seqmat, bs, rows, bs, cols, vals, INSERT_VALUES));
55c4762a1bSJed Brown     }
56c4762a1bSJed Brown     /* off-diagonal blocks */
57c4762a1bSJed Brown     for (i = 0; i < bs * bs; i++) vals[i] = -1.0;
58c4762a1bSJed Brown     for (i = 0; i < (n - 1) * bs; i += bs) {
599371c9d4SSatish Balay       for (j = 0; j < bs; j++) {
609371c9d4SSatish Balay         rows[j] = i + j;
619371c9d4SSatish Balay         cols[j] = i + bs + j;
629371c9d4SSatish Balay       }
639566063dSJacob Faibussowitsch       PetscCall(MatSetValues(seqmat, bs, rows, bs, cols, vals, INSERT_VALUES));
64c4762a1bSJed Brown     }
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch     PetscCall(PetscFree3(vals, rows, cols));
67c4762a1bSJed Brown   }
689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(seqmat, MAT_FINAL_ASSEMBLY));
699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(seqmat, MAT_FINAL_ASSEMBLY));
70dd400576SPatrick Sanan   if (rank == 0) {
719566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] seqmat:\n", rank));
729566063dSJacob Faibussowitsch     PetscCall(MatView(seqmat, PETSC_VIEWER_STDOUT_SELF));
73c4762a1bSJed Brown   }
74c4762a1bSJed Brown 
759566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD, seqmat, PETSC_DECIDE, MAT_INITIAL_MATRIX, &mpimat));
769566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD, seqmat, PETSC_DECIDE, MAT_REUSE_MATRIX, &mpimat));
779566063dSJacob Faibussowitsch   PetscCall(MatView(mpimat, PETSC_VIEWER_STDOUT_WORLD));
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&seqmat));
809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mpimat));
819566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
82b122ec5aSJacob Faibussowitsch   return 0;
83c4762a1bSJed Brown }
84c4762a1bSJed Brown 
85c4762a1bSJed Brown /*TEST
86c4762a1bSJed Brown 
87c4762a1bSJed Brown    test:
88c4762a1bSJed Brown       nsize: 3
89c4762a1bSJed Brown 
90c4762a1bSJed Brown    test:
91c4762a1bSJed Brown       suffix: 2
92c4762a1bSJed Brown       nsize: 3
93c4762a1bSJed Brown       args: -mat_type baij
94c4762a1bSJed Brown 
95c4762a1bSJed Brown    test:
96c4762a1bSJed Brown       suffix: 3
97c4762a1bSJed Brown       nsize: 3
98c4762a1bSJed Brown       args: -mat_type baij -bs 2
99c4762a1bSJed Brown 
100c4762a1bSJed Brown TEST*/
101