1c4762a1bSJed Brown static char help[] = "Tests MPI parallel matrix creation. Test MatCreateRedundantMatrix() \n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown Mat C, Credundant;
8c4762a1bSJed Brown MatInfo info;
9c4762a1bSJed Brown PetscMPIInt rank, size, subsize;
10c4762a1bSJed Brown PetscInt i, j, m = 3, n = 2, low, high, iglobal;
11c4762a1bSJed Brown PetscInt Ii, J, ldim, nsubcomms;
12c4762a1bSJed Brown PetscBool flg_info, flg_mat;
13c4762a1bSJed Brown PetscScalar v, one = 1.0;
14c4762a1bSJed Brown Vec x, y;
15c4762a1bSJed Brown MPI_Comm subcomm;
16c4762a1bSJed Brown
17327415f7SBarry Smith PetscFunctionBeginUser;
18*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
22c4762a1bSJed Brown n = 2 * size;
23c4762a1bSJed Brown
249566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
259566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
269566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C));
279566063dSJacob Faibussowitsch PetscCall(MatSetUp(C));
28c4762a1bSJed Brown
29c4762a1bSJed Brown /* Create the matrix for the five point stencil, YET AGAIN */
30c4762a1bSJed Brown for (i = 0; i < m; i++) {
31c4762a1bSJed Brown for (j = 2 * rank; j < 2 * rank + 2; j++) {
329371c9d4SSatish Balay v = -1.0;
339371c9d4SSatish Balay Ii = j + n * i;
349371c9d4SSatish Balay if (i > 0) {
359371c9d4SSatish Balay J = Ii - n;
369371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
379371c9d4SSatish Balay }
389371c9d4SSatish Balay if (i < m - 1) {
399371c9d4SSatish Balay J = Ii + n;
409371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
419371c9d4SSatish Balay }
429371c9d4SSatish Balay if (j > 0) {
439371c9d4SSatish Balay J = Ii - 1;
449371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
459371c9d4SSatish Balay }
469371c9d4SSatish Balay if (j < n - 1) {
479371c9d4SSatish Balay J = Ii + 1;
489371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
499371c9d4SSatish Balay }
509371c9d4SSatish Balay v = 4.0;
519371c9d4SSatish Balay PetscCall(MatSetValues(C, 1, &Ii, 1, &Ii, &v, INSERT_VALUES));
52c4762a1bSJed Brown }
53c4762a1bSJed Brown }
54c4762a1bSJed Brown
55c4762a1bSJed Brown /* Add extra elements (to illustrate variants of MatGetInfo) */
569371c9d4SSatish Balay Ii = n;
579371c9d4SSatish Balay J = n - 2;
589371c9d4SSatish Balay v = 100.0;
599566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
609371c9d4SSatish Balay Ii = n - 2;
619371c9d4SSatish Balay J = n;
629371c9d4SSatish Balay v = 100.0;
639566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &Ii, 1, &J, &v, INSERT_VALUES));
64c4762a1bSJed Brown
659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
67c4762a1bSJed Brown
68c4762a1bSJed Brown /* Form vectors */
699566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(C, &x, &y));
709566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &ldim));
719566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(x, &low, &high));
72c4762a1bSJed Brown for (i = 0; i < ldim; i++) {
73c4762a1bSJed Brown iglobal = i + low;
74c4762a1bSJed Brown v = one * ((PetscReal)i) + 100.0 * rank;
759566063dSJacob Faibussowitsch PetscCall(VecSetValues(x, 1, &iglobal, &v, INSERT_VALUES));
76c4762a1bSJed Brown }
779566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x));
789566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x));
79c4762a1bSJed Brown
809566063dSJacob Faibussowitsch PetscCall(MatMult(C, x, y));
81c4762a1bSJed Brown
829566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-view_info", &flg_info));
83c4762a1bSJed Brown if (flg_info) {
849566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO));
859566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
86c4762a1bSJed Brown
879566063dSJacob Faibussowitsch PetscCall(MatGetInfo(C, MAT_GLOBAL_SUM, &info));
889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "matrix information (global sums):\nnonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n", (PetscInt)info.nz_used, (PetscInt)info.nz_allocated));
899566063dSJacob Faibussowitsch PetscCall(MatGetInfo(C, MAT_GLOBAL_MAX, &info));
909566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "matrix information (global max):\nnonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n", (PetscInt)info.nz_used, (PetscInt)info.nz_allocated));
91c4762a1bSJed Brown }
92c4762a1bSJed Brown
939566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-view_mat", &flg_mat));
941baa6e33SBarry Smith if (flg_mat) PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
95c4762a1bSJed Brown
96c4762a1bSJed Brown /* Test MatCreateRedundantMatrix() */
97c4762a1bSJed Brown nsubcomms = size;
989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nsubcomms", &nsubcomms, NULL));
999566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(C, nsubcomms, MPI_COMM_NULL, MAT_INITIAL_MATRIX, &Credundant));
1009566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(C, nsubcomms, MPI_COMM_NULL, MAT_REUSE_MATRIX, &Credundant));
101c4762a1bSJed Brown
1029566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)Credundant, &subcomm));
1039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(subcomm, &subsize));
104c4762a1bSJed Brown
105c4762a1bSJed Brown if (subsize == 2 && flg_mat) {
1069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_(subcomm), "\n[%d] Credundant:\n", rank));
1079566063dSJacob Faibussowitsch PetscCall(MatView(Credundant, PETSC_VIEWER_STDOUT_(subcomm)));
108c4762a1bSJed Brown }
1099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Credundant));
110c4762a1bSJed Brown
111c4762a1bSJed Brown /* Test MatCreateRedundantMatrix() with user-provided subcomm */
112c4762a1bSJed Brown {
113c4762a1bSJed Brown PetscSubcomm psubcomm;
114c4762a1bSJed Brown
1159566063dSJacob Faibussowitsch PetscCall(PetscSubcommCreate(PETSC_COMM_WORLD, &psubcomm));
1169566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetNumber(psubcomm, nsubcomms));
1179566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetType(psubcomm, PETSC_SUBCOMM_CONTIGUOUS));
118c4762a1bSJed Brown /* enable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */
1199566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetFromOptions(psubcomm));
120c4762a1bSJed Brown
1219566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(C, nsubcomms, PetscSubcommChild(psubcomm), MAT_INITIAL_MATRIX, &Credundant));
1229566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(C, nsubcomms, PetscSubcommChild(psubcomm), MAT_REUSE_MATRIX, &Credundant));
123c4762a1bSJed Brown
1249566063dSJacob Faibussowitsch PetscCall(PetscSubcommDestroy(&psubcomm));
1259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Credundant));
126c4762a1bSJed Brown }
127c4762a1bSJed Brown
1289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y));
1309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
1319566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
132b122ec5aSJacob Faibussowitsch return 0;
133c4762a1bSJed Brown }
134c4762a1bSJed Brown
135c4762a1bSJed Brown /*TEST
136c4762a1bSJed Brown
137c4762a1bSJed Brown test:
138c4762a1bSJed Brown nsize: 3
139c4762a1bSJed Brown args: -view_info
140c4762a1bSJed Brown
141c4762a1bSJed Brown test:
142c4762a1bSJed Brown suffix: 2
143c4762a1bSJed Brown nsize: 3
144c4762a1bSJed Brown args: -nsubcomms 2 -view_mat -psubcomm_type interlaced
145c4762a1bSJed Brown
146c4762a1bSJed Brown test:
147c4762a1bSJed Brown suffix: 3
148c4762a1bSJed Brown nsize: 3
149c4762a1bSJed Brown args: -nsubcomms 2 -view_mat -psubcomm_type contiguous
150c4762a1bSJed Brown
151c4762a1bSJed Brown test:
152c4762a1bSJed Brown suffix: 3_baij
153c4762a1bSJed Brown nsize: 3
154c4762a1bSJed Brown args: -mat_type baij -nsubcomms 2 -view_mat
155c4762a1bSJed Brown
156c4762a1bSJed Brown test:
157c4762a1bSJed Brown suffix: 3_sbaij
158c4762a1bSJed Brown nsize: 3
159c4762a1bSJed Brown args: -mat_type sbaij -nsubcomms 2 -view_mat
160c4762a1bSJed Brown
161c4762a1bSJed Brown test:
162c4762a1bSJed Brown suffix: 3_dense
163c4762a1bSJed Brown nsize: 3
164c4762a1bSJed Brown args: -mat_type dense -nsubcomms 2 -view_mat
165c4762a1bSJed Brown
166c4762a1bSJed Brown test:
167c4762a1bSJed Brown suffix: 4_baij
168c4762a1bSJed Brown nsize: 3
169c4762a1bSJed Brown args: -mat_type baij -nsubcomms 2 -view_mat -psubcomm_type interlaced
170c4762a1bSJed Brown
171c4762a1bSJed Brown test:
172c4762a1bSJed Brown suffix: 4_sbaij
173c4762a1bSJed Brown nsize: 3
174c4762a1bSJed Brown args: -mat_type sbaij -nsubcomms 2 -view_mat -psubcomm_type interlaced
175c4762a1bSJed Brown
176c4762a1bSJed Brown test:
177c4762a1bSJed Brown suffix: 4_dense
178c4762a1bSJed Brown nsize: 3
179c4762a1bSJed Brown args: -mat_type dense -nsubcomms 2 -view_mat -psubcomm_type interlaced
180c4762a1bSJed Brown
181c4762a1bSJed Brown TEST*/
182