xref: /petsc/src/mat/tests/ex86.c (revision ffeef943c8ee50edff320d8a3135bb0c94853e4c)
1 static char help[] = "Testing MatCreateMPIMatConcatenateSeqMat().\n\n";
2 
3 #include <petscmat.h>
4 int main(int argc, char **argv)
5 {
6   Mat         seqmat, mpimat;
7   PetscMPIInt rank;
8   PetscScalar value[3], *vals;
9   PetscInt    i, col[3], n = 5, bs = 1;
10 
11   PetscFunctionBeginUser;
12   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
13   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
14   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
15 
16   /* Create seqaij matrices of size (n+rank) by n */
17   PetscCall(MatCreate(PETSC_COMM_SELF, &seqmat));
18   PetscCall(MatSetSizes(seqmat, (n + rank) * bs, PETSC_DECIDE, PETSC_DECIDE, n * bs));
19   PetscCall(MatSetFromOptions(seqmat));
20   PetscCall(MatSeqAIJSetPreallocation(seqmat, 3 * bs, NULL));
21   PetscCall(MatSeqBAIJSetPreallocation(seqmat, bs, 3, NULL));
22 
23   if (bs == 1) {
24     value[0] = -1.0;
25     value[1] = 2.0;
26     value[2] = -1.0;
27     for (i = 1; i < n - 1; i++) {
28       col[0] = i - 1;
29       col[1] = i;
30       col[2] = i + 1;
31       PetscCall(MatSetValues(seqmat, 1, &i, 3, col, value, INSERT_VALUES));
32     }
33     i      = n - 1;
34     col[0] = n - 2;
35     col[1] = n - 1;
36     PetscCall(MatSetValues(seqmat, 1, &i, 2, col, value, INSERT_VALUES));
37 
38     i        = 0;
39     col[0]   = 0;
40     col[1]   = 1;
41     value[0] = 2.0;
42     value[1] = -1.0;
43     PetscCall(MatSetValues(seqmat, 1, &i, 2, col, value, INSERT_VALUES));
44   } else {
45     PetscInt *rows, *cols, j;
46     PetscCall(PetscMalloc3(bs * bs, &vals, bs, &rows, bs, &cols));
47     /* diagonal blocks */
48     for (i = 0; i < bs * bs; i++) vals[i] = 2.0;
49     for (i = 0; i < n * bs; i += bs) {
50       for (j = 0; j < bs; j++) {
51         rows[j] = i + j;
52         cols[j] = i + j;
53       }
54       PetscCall(MatSetValues(seqmat, bs, rows, bs, cols, vals, INSERT_VALUES));
55     }
56     /* off-diagonal blocks */
57     for (i = 0; i < bs * bs; i++) vals[i] = -1.0;
58     for (i = 0; i < (n - 1) * bs; i += bs) {
59       for (j = 0; j < bs; j++) {
60         rows[j] = i + j;
61         cols[j] = i + bs + j;
62       }
63       PetscCall(MatSetValues(seqmat, bs, rows, bs, cols, vals, INSERT_VALUES));
64     }
65 
66     PetscCall(PetscFree3(vals, rows, cols));
67   }
68   PetscCall(MatAssemblyBegin(seqmat, MAT_FINAL_ASSEMBLY));
69   PetscCall(MatAssemblyEnd(seqmat, MAT_FINAL_ASSEMBLY));
70   if (rank == 0) {
71     PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] seqmat:\n", rank));
72     PetscCall(MatView(seqmat, PETSC_VIEWER_STDOUT_SELF));
73   }
74 
75   PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD, seqmat, PETSC_DECIDE, MAT_INITIAL_MATRIX, &mpimat));
76   PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD, seqmat, PETSC_DECIDE, MAT_REUSE_MATRIX, &mpimat));
77   PetscCall(MatView(mpimat, PETSC_VIEWER_STDOUT_WORLD));
78 
79   PetscCall(MatDestroy(&seqmat));
80   PetscCall(MatDestroy(&mpimat));
81   PetscCall(PetscFinalize());
82   return 0;
83 }
84 
85 /*TEST
86 
87    test:
88       nsize: 3
89 
90    test:
91       suffix: 2
92       nsize: 3
93       args: -mat_type baij
94 
95    test:
96       suffix: 3
97       nsize: 3
98       args: -mat_type baij -bs 2
99 
100 TEST*/
101