xref: /petsc/src/mat/tests/ex86.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Testing MatCreateMPIMatConcatenateSeqMat().\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4*9371c9d4SSatish Balay int main(int argc, char **argv) {
5c4762a1bSJed Brown   Mat         seqmat, mpimat;
6c4762a1bSJed Brown   PetscMPIInt rank;
7c4762a1bSJed Brown   PetscScalar value[3], *vals;
8c4762a1bSJed Brown   PetscInt    i, col[3], n = 5, bs = 1;
9c4762a1bSJed Brown 
10327415f7SBarry Smith   PetscFunctionBeginUser;
119566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   /* Create seqaij matrices of size (n+rank) by n */
169566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &seqmat));
179566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(seqmat, (n + rank) * bs, PETSC_DECIDE, PETSC_DECIDE, n * bs));
189566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(seqmat));
199566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(seqmat, 3 * bs, NULL));
209566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(seqmat, bs, 3, NULL));
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   if (bs == 1) {
23*9371c9d4SSatish Balay     value[0] = -1.0;
24*9371c9d4SSatish Balay     value[1] = 2.0;
25*9371c9d4SSatish Balay     value[2] = -1.0;
26c4762a1bSJed Brown     for (i = 1; i < n - 1; i++) {
27*9371c9d4SSatish Balay       col[0] = i - 1;
28*9371c9d4SSatish Balay       col[1] = i;
29*9371c9d4SSatish Balay       col[2] = i + 1;
309566063dSJacob Faibussowitsch       PetscCall(MatSetValues(seqmat, 1, &i, 3, col, value, INSERT_VALUES));
31c4762a1bSJed Brown     }
32*9371c9d4SSatish Balay     i      = n - 1;
33*9371c9d4SSatish Balay     col[0] = n - 2;
34*9371c9d4SSatish Balay     col[1] = n - 1;
359566063dSJacob Faibussowitsch     PetscCall(MatSetValues(seqmat, 1, &i, 2, col, value, INSERT_VALUES));
36c4762a1bSJed Brown 
37*9371c9d4SSatish Balay     i        = 0;
38*9371c9d4SSatish Balay     col[0]   = 0;
39*9371c9d4SSatish Balay     col[1]   = 1;
40*9371c9d4SSatish Balay     value[0] = 2.0;
41*9371c9d4SSatish Balay     value[1] = -1.0;
429566063dSJacob Faibussowitsch     PetscCall(MatSetValues(seqmat, 1, &i, 2, col, value, INSERT_VALUES));
43c4762a1bSJed Brown   } else {
44c4762a1bSJed Brown     PetscInt *rows, *cols, j;
459566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(bs * bs, &vals, bs, &rows, bs, &cols));
46c4762a1bSJed Brown     /* diagonal blocks */
47c4762a1bSJed Brown     for (i = 0; i < bs * bs; i++) vals[i] = 2.0;
48c4762a1bSJed Brown     for (i = 0; i < n * bs; i += bs) {
49*9371c9d4SSatish Balay       for (j = 0; j < bs; j++) {
50*9371c9d4SSatish Balay         rows[j] = i + j;
51*9371c9d4SSatish Balay         cols[j] = i + j;
52*9371c9d4SSatish Balay       }
539566063dSJacob Faibussowitsch       PetscCall(MatSetValues(seqmat, bs, rows, bs, cols, vals, INSERT_VALUES));
54c4762a1bSJed Brown     }
55c4762a1bSJed Brown     /* off-diagonal blocks */
56c4762a1bSJed Brown     for (i = 0; i < bs * bs; i++) vals[i] = -1.0;
57c4762a1bSJed Brown     for (i = 0; i < (n - 1) * bs; i += bs) {
58*9371c9d4SSatish Balay       for (j = 0; j < bs; j++) {
59*9371c9d4SSatish Balay         rows[j] = i + j;
60*9371c9d4SSatish Balay         cols[j] = i + bs + j;
61*9371c9d4SSatish Balay       }
629566063dSJacob Faibussowitsch       PetscCall(MatSetValues(seqmat, bs, rows, bs, cols, vals, INSERT_VALUES));
63c4762a1bSJed Brown     }
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch     PetscCall(PetscFree3(vals, rows, cols));
66c4762a1bSJed Brown   }
679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(seqmat, MAT_FINAL_ASSEMBLY));
689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(seqmat, MAT_FINAL_ASSEMBLY));
69dd400576SPatrick Sanan   if (rank == 0) {
709566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] seqmat:\n", rank));
719566063dSJacob Faibussowitsch     PetscCall(MatView(seqmat, PETSC_VIEWER_STDOUT_SELF));
72c4762a1bSJed Brown   }
73c4762a1bSJed Brown 
749566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD, seqmat, PETSC_DECIDE, MAT_INITIAL_MATRIX, &mpimat));
759566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD, seqmat, PETSC_DECIDE, MAT_REUSE_MATRIX, &mpimat));
769566063dSJacob Faibussowitsch   PetscCall(MatView(mpimat, PETSC_VIEWER_STDOUT_WORLD));
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&seqmat));
799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mpimat));
809566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
81b122ec5aSJacob Faibussowitsch   return 0;
82c4762a1bSJed Brown }
83c4762a1bSJed Brown 
84c4762a1bSJed Brown /*TEST
85c4762a1bSJed Brown 
86c4762a1bSJed Brown    test:
87c4762a1bSJed Brown       nsize: 3
88c4762a1bSJed Brown 
89c4762a1bSJed Brown    test:
90c4762a1bSJed Brown       suffix: 2
91c4762a1bSJed Brown       nsize: 3
92c4762a1bSJed Brown       args: -mat_type baij
93c4762a1bSJed Brown 
94c4762a1bSJed Brown    test:
95c4762a1bSJed Brown       suffix: 3
96c4762a1bSJed Brown       nsize: 3
97c4762a1bSJed Brown       args: -mat_type baij -bs 2
98c4762a1bSJed Brown 
99c4762a1bSJed Brown TEST*/
100