xref: /petsc/src/ksp/pc/impls/telescope/telescope_dmda.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
11e07b27eSBarry Smith 
2120bdd93SDave May #include <petsc/private/matimpl.h>
31e07b27eSBarry Smith #include <petsc/private/pcimpl.h>
45e897e82SDave May #include <petsc/private/dmimpl.h>
51e07b27eSBarry Smith #include <petscksp.h> /*I "petscksp.h" I*/
61e07b27eSBarry Smith #include <petscdm.h>
71e07b27eSBarry Smith #include <petscdmda.h>
81e07b27eSBarry Smith 
9575a0592SBarry Smith #include "../src/ksp/pc/impls/telescope/telescope.h"
101e07b27eSBarry Smith 
11bf00f589SPatrick Sanan static PetscBool  cited      = PETSC_FALSE;
129371c9d4SSatish Balay static const char citation[] = "@inproceedings{MaySananRuppKnepleySmith2016,\n"
13bf00f589SPatrick Sanan                                "  title     = {Extreme-Scale Multigrid Components within PETSc},\n"
14bf00f589SPatrick Sanan                                "  author    = {Dave A. May and Patrick Sanan and Karl Rupp and Matthew G. Knepley and Barry F. Smith},\n"
15bf00f589SPatrick Sanan                                "  booktitle = {Proceedings of the Platform for Advanced Scientific Computing Conference},\n"
16bf00f589SPatrick Sanan                                "  series    = {PASC '16},\n"
17bf00f589SPatrick Sanan                                "  isbn      = {978-1-4503-4126-4},\n"
18bf00f589SPatrick Sanan                                "  location  = {Lausanne, Switzerland},\n"
19bf00f589SPatrick Sanan                                "  pages     = {5:1--5:12},\n"
20bf00f589SPatrick Sanan                                "  articleno = {5},\n"
21bf00f589SPatrick Sanan                                "  numpages  = {12},\n"
22a8d69d7bSBarry Smith                                "  url       = {https://doi.acm.org/10.1145/2929908.2929913},\n"
23bf00f589SPatrick Sanan                                "  doi       = {10.1145/2929908.2929913},\n"
24bf00f589SPatrick Sanan                                "  acmid     = {2929913},\n"
25bf00f589SPatrick Sanan                                "  publisher = {ACM},\n"
26bf00f589SPatrick Sanan                                "  address   = {New York, NY, USA},\n"
27bf00f589SPatrick Sanan                                "  keywords  = {GPU, HPC, agglomeration, coarse-level solver, multigrid, parallel computing, preconditioning},\n"
28bf00f589SPatrick Sanan                                "  year      = {2016}\n"
29bf00f589SPatrick Sanan                                "}\n";
30bf00f589SPatrick Sanan 
319371c9d4SSatish Balay static PetscErrorCode _DMDADetermineRankFromGlobalIJK(PetscInt dim, PetscInt i, PetscInt j, PetscInt k, PetscInt Mp, PetscInt Np, PetscInt Pp, PetscInt start_i[], PetscInt start_j[], PetscInt start_k[], PetscInt span_i[], PetscInt span_j[], PetscInt span_k[], PetscMPIInt *_pi, PetscMPIInt *_pj, PetscMPIInt *_pk, PetscMPIInt *rank_re) {
321e07b27eSBarry Smith   PetscInt pi, pj, pk, n;
331e07b27eSBarry Smith 
341e07b27eSBarry Smith   PetscFunctionBegin;
35137d0469SJed Brown   *rank_re = -1;
36137d0469SJed Brown   if (_pi) *_pi = -1;
37137d0469SJed Brown   if (_pj) *_pj = -1;
38137d0469SJed Brown   if (_pk) *_pk = -1;
391e07b27eSBarry Smith   pi = pj = pk = -1;
401e07b27eSBarry Smith   if (_pi) {
411e07b27eSBarry Smith     for (n = 0; n < Mp; n++) {
421e07b27eSBarry Smith       if ((i >= start_i[n]) && (i < start_i[n] + span_i[n])) {
431e07b27eSBarry Smith         pi = n;
441e07b27eSBarry Smith         break;
451e07b27eSBarry Smith       }
461e07b27eSBarry Smith     }
4763a3b9bcSJacob Faibussowitsch     PetscCheck(pi != -1, PETSC_COMM_SELF, PETSC_ERR_USER, "[dmda-ijk] pi cannot be determined : range %" PetscInt_FMT ", val %" PetscInt_FMT, Mp, i);
481e07b27eSBarry Smith     *_pi = pi;
491e07b27eSBarry Smith   }
501e07b27eSBarry Smith 
511e07b27eSBarry Smith   if (_pj) {
521e07b27eSBarry Smith     for (n = 0; n < Np; n++) {
531e07b27eSBarry Smith       if ((j >= start_j[n]) && (j < start_j[n] + span_j[n])) {
541e07b27eSBarry Smith         pj = n;
551e07b27eSBarry Smith         break;
561e07b27eSBarry Smith       }
571e07b27eSBarry Smith     }
5863a3b9bcSJacob Faibussowitsch     PetscCheck(pj != -1, PETSC_COMM_SELF, PETSC_ERR_USER, "[dmda-ijk] pj cannot be determined : range %" PetscInt_FMT ", val %" PetscInt_FMT, Np, j);
591e07b27eSBarry Smith     *_pj = pj;
601e07b27eSBarry Smith   }
611e07b27eSBarry Smith 
621e07b27eSBarry Smith   if (_pk) {
631e07b27eSBarry Smith     for (n = 0; n < Pp; n++) {
641e07b27eSBarry Smith       if ((k >= start_k[n]) && (k < start_k[n] + span_k[n])) {
651e07b27eSBarry Smith         pk = n;
661e07b27eSBarry Smith         break;
671e07b27eSBarry Smith       }
681e07b27eSBarry Smith     }
6963a3b9bcSJacob Faibussowitsch     PetscCheck(pk != -1, PETSC_COMM_SELF, PETSC_ERR_USER, "[dmda-ijk] pk cannot be determined : range %" PetscInt_FMT ", val %" PetscInt_FMT, Pp, k);
701e07b27eSBarry Smith     *_pk = pk;
711e07b27eSBarry Smith   }
721e07b27eSBarry Smith 
731e07b27eSBarry Smith   switch (dim) {
749371c9d4SSatish Balay   case 1: *rank_re = pi; break;
759371c9d4SSatish Balay   case 2: *rank_re = pi + pj * Mp; break;
769371c9d4SSatish Balay   case 3: *rank_re = pi + pj * Mp + pk * (Mp * Np); break;
771e07b27eSBarry Smith   }
781e07b27eSBarry Smith   PetscFunctionReturn(0);
791e07b27eSBarry Smith }
801e07b27eSBarry Smith 
819371c9d4SSatish Balay static PetscErrorCode _DMDADetermineGlobalS0(PetscInt dim, PetscMPIInt rank_re, PetscInt Mp_re, PetscInt Np_re, PetscInt Pp_re, PetscInt range_i_re[], PetscInt range_j_re[], PetscInt range_k_re[], PetscInt *s0) {
82c6a0d831SBarry Smith   PetscInt i, j, k, start_IJK = 0;
831e07b27eSBarry Smith   PetscInt rank_ijk;
841e07b27eSBarry Smith 
851e07b27eSBarry Smith   PetscFunctionBegin;
861e07b27eSBarry Smith   switch (dim) {
871e07b27eSBarry Smith   case 1:
881e07b27eSBarry Smith     for (i = 0; i < Mp_re; i++) {
891e07b27eSBarry Smith       rank_ijk = i;
909371c9d4SSatish Balay       if (rank_ijk < rank_re) { start_IJK += range_i_re[i]; }
911e07b27eSBarry Smith     }
921e07b27eSBarry Smith     break;
931e07b27eSBarry Smith   case 2:
941e07b27eSBarry Smith     for (j = 0; j < Np_re; j++) {
951e07b27eSBarry Smith       for (i = 0; i < Mp_re; i++) {
961e07b27eSBarry Smith         rank_ijk = i + j * Mp_re;
979371c9d4SSatish Balay         if (rank_ijk < rank_re) { start_IJK += range_i_re[i] * range_j_re[j]; }
981e07b27eSBarry Smith       }
991e07b27eSBarry Smith     }
1001e07b27eSBarry Smith     break;
1011e07b27eSBarry Smith   case 3:
1021e07b27eSBarry Smith     for (k = 0; k < Pp_re; k++) {
1031e07b27eSBarry Smith       for (j = 0; j < Np_re; j++) {
1041e07b27eSBarry Smith         for (i = 0; i < Mp_re; i++) {
1051e07b27eSBarry Smith           rank_ijk = i + j * Mp_re + k * Mp_re * Np_re;
1069371c9d4SSatish Balay           if (rank_ijk < rank_re) { start_IJK += range_i_re[i] * range_j_re[j] * range_k_re[k]; }
1071e07b27eSBarry Smith         }
1081e07b27eSBarry Smith       }
1091e07b27eSBarry Smith     }
1101e07b27eSBarry Smith     break;
1111e07b27eSBarry Smith   }
1121e07b27eSBarry Smith   *s0 = start_IJK;
1131e07b27eSBarry Smith   PetscFunctionReturn(0);
1141e07b27eSBarry Smith }
1151e07b27eSBarry Smith 
1169371c9d4SSatish Balay static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors2d(PC_Telescope sred, DM dm, DM subdm) {
1171e07b27eSBarry Smith   DM         cdm;
1181e07b27eSBarry Smith   Vec        coor, coor_natural, perm_coors;
1191e07b27eSBarry Smith   PetscInt   i, j, si, sj, ni, nj, M, N, Ml, Nl, c, nidx;
1201e07b27eSBarry Smith   PetscInt  *fine_indices;
1211e07b27eSBarry Smith   IS         is_fine, is_local;
1221e07b27eSBarry Smith   VecScatter sctx;
1231e07b27eSBarry Smith 
1241e07b27eSBarry Smith   PetscFunctionBegin;
1259566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(dm, &coor));
1261e07b27eSBarry Smith   if (!coor) return (0);
127*48a46eb9SPierre Jolivet   if (PCTelescope_isActiveRank(sred)) PetscCall(DMDASetUniformCoordinates(subdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
1281e07b27eSBarry Smith   /* Get the coordinate vector from the distributed array */
1299566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
1309566063dSJacob Faibussowitsch   PetscCall(DMDACreateNaturalVector(cdm, &coor_natural));
1311e07b27eSBarry Smith 
1329566063dSJacob Faibussowitsch   PetscCall(DMDAGlobalToNaturalBegin(cdm, coor, INSERT_VALUES, coor_natural));
1339566063dSJacob Faibussowitsch   PetscCall(DMDAGlobalToNaturalEnd(cdm, coor, INSERT_VALUES, coor_natural));
1341e07b27eSBarry Smith 
1351e07b27eSBarry Smith   /* get indices of the guys I want to grab */
1369566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
13757f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
1389566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(subdm, &si, &sj, NULL, &ni, &nj, NULL));
13915dd08bcSBarry Smith     Ml = ni;
14015dd08bcSBarry Smith     Nl = nj;
1411e07b27eSBarry Smith   } else {
142c41e779fSDave May     si = sj = 0;
143c41e779fSDave May     ni = nj = 0;
1443ac26c5eSBarry Smith     Ml = Nl = 0;
1451e07b27eSBarry Smith   }
1461e07b27eSBarry Smith 
1479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Ml * Nl * 2, &fine_indices));
1481e07b27eSBarry Smith   c = 0;
14957f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
1501e07b27eSBarry Smith     for (j = sj; j < sj + nj; j++) {
1511e07b27eSBarry Smith       for (i = si; i < si + ni; i++) {
1521e07b27eSBarry Smith         nidx                = (i) + (j)*M;
1531e07b27eSBarry Smith         fine_indices[c]     = 2 * nidx;
1541e07b27eSBarry Smith         fine_indices[c + 1] = 2 * nidx + 1;
1551e07b27eSBarry Smith         c                   = c + 2;
1561e07b27eSBarry Smith       }
1571e07b27eSBarry Smith     }
15863a3b9bcSJacob Faibussowitsch     PetscCheck(c == Ml * Nl * 2, PETSC_COMM_SELF, PETSC_ERR_PLIB, "c %" PetscInt_FMT " should equal 2 * Ml %" PetscInt_FMT " * Nl %" PetscInt_FMT, c, Ml, Nl);
1591e07b27eSBarry Smith   }
1601e07b27eSBarry Smith 
1611e07b27eSBarry Smith   /* generate scatter */
1629566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), Ml * Nl * 2, fine_indices, PETSC_USE_POINTER, &is_fine));
1639566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, Ml * Nl * 2, 0, 1, &is_local));
1641e07b27eSBarry Smith 
1651e07b27eSBarry Smith   /* scatter */
1669566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &perm_coors));
1679566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(perm_coors, PETSC_DECIDE, Ml * Nl * 2));
1689566063dSJacob Faibussowitsch   PetscCall(VecSetType(perm_coors, VECSEQ));
1691e07b27eSBarry Smith 
1709566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(coor_natural, is_fine, perm_coors, is_local, &sctx));
1719566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sctx, coor_natural, perm_coors, INSERT_VALUES, SCATTER_FORWARD));
1729566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sctx, coor_natural, perm_coors, INSERT_VALUES, SCATTER_FORWARD));
1731e07b27eSBarry Smith   /* access */
17457f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
1751e07b27eSBarry Smith     Vec                _coors;
1761e07b27eSBarry Smith     const PetscScalar *LA_perm;
1771e07b27eSBarry Smith     PetscScalar       *LA_coors;
1781e07b27eSBarry Smith 
1799566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(subdm, &_coors));
1809566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(perm_coors, &LA_perm));
1819566063dSJacob Faibussowitsch     PetscCall(VecGetArray(_coors, &LA_coors));
1829371c9d4SSatish Balay     for (i = 0; i < Ml * Nl * 2; i++) { LA_coors[i] = LA_perm[i]; }
1839566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(_coors, &LA_coors));
1849566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(perm_coors, &LA_perm));
1851e07b27eSBarry Smith   }
1861e07b27eSBarry Smith 
1871e07b27eSBarry Smith   /* update local coords */
18857f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
1891e07b27eSBarry Smith     DM  _dmc;
1901e07b27eSBarry Smith     Vec _coors, _coors_local;
1919566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(subdm, &_dmc));
1929566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(subdm, &_coors));
1939566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(subdm, &_coors_local));
1949566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(_dmc, _coors, INSERT_VALUES, _coors_local));
1959566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(_dmc, _coors, INSERT_VALUES, _coors_local));
1961e07b27eSBarry Smith   }
1979566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&sctx));
1989566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_fine));
1999566063dSJacob Faibussowitsch   PetscCall(PetscFree(fine_indices));
2009566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_local));
2019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&perm_coors));
2029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coor_natural));
2031e07b27eSBarry Smith   PetscFunctionReturn(0);
2041e07b27eSBarry Smith }
2051e07b27eSBarry Smith 
2069371c9d4SSatish Balay static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors3d(PC_Telescope sred, DM dm, DM subdm) {
2071e07b27eSBarry Smith   DM         cdm;
2081e07b27eSBarry Smith   Vec        coor, coor_natural, perm_coors;
2091e07b27eSBarry Smith   PetscInt   i, j, k, si, sj, sk, ni, nj, nk, M, N, P, Ml, Nl, Pl, c, nidx;
2101e07b27eSBarry Smith   PetscInt  *fine_indices;
2111e07b27eSBarry Smith   IS         is_fine, is_local;
2121e07b27eSBarry Smith   VecScatter sctx;
2131e07b27eSBarry Smith 
2141e07b27eSBarry Smith   PetscFunctionBegin;
2159566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(dm, &coor));
2161e07b27eSBarry Smith   if (!coor) PetscFunctionReturn(0);
2171e07b27eSBarry Smith 
218*48a46eb9SPierre Jolivet   if (PCTelescope_isActiveRank(sred)) PetscCall(DMDASetUniformCoordinates(subdm, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
2191e07b27eSBarry Smith 
2201e07b27eSBarry Smith   /* Get the coordinate vector from the distributed array */
2219566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
2229566063dSJacob Faibussowitsch   PetscCall(DMDACreateNaturalVector(cdm, &coor_natural));
2239566063dSJacob Faibussowitsch   PetscCall(DMDAGlobalToNaturalBegin(cdm, coor, INSERT_VALUES, coor_natural));
2249566063dSJacob Faibussowitsch   PetscCall(DMDAGlobalToNaturalEnd(cdm, coor, INSERT_VALUES, coor_natural));
2251e07b27eSBarry Smith 
2261e07b27eSBarry Smith   /* get indices of the guys I want to grab */
2279566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, NULL, &M, &N, &P, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
2281e07b27eSBarry Smith 
22957f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
2309566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(subdm, &si, &sj, &sk, &ni, &nj, &nk));
231553d0ae9SBarry Smith     Ml = ni;
232553d0ae9SBarry Smith     Nl = nj;
233553d0ae9SBarry Smith     Pl = nk;
2341e07b27eSBarry Smith   } else {
235c41e779fSDave May     si = sj = sk = 0;
236c41e779fSDave May     ni = nj = nk = 0;
2373ac26c5eSBarry Smith     Ml = Nl = Pl = 0;
2381e07b27eSBarry Smith   }
2391e07b27eSBarry Smith 
2409566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Ml * Nl * Pl * 3, &fine_indices));
2411e07b27eSBarry Smith 
2421e07b27eSBarry Smith   c = 0;
24357f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
2441e07b27eSBarry Smith     for (k = sk; k < sk + nk; k++) {
2451e07b27eSBarry Smith       for (j = sj; j < sj + nj; j++) {
2461e07b27eSBarry Smith         for (i = si; i < si + ni; i++) {
2471e07b27eSBarry Smith           nidx                = (i) + (j)*M + (k)*M * N;
2481e07b27eSBarry Smith           fine_indices[c]     = 3 * nidx;
2491e07b27eSBarry Smith           fine_indices[c + 1] = 3 * nidx + 1;
2501e07b27eSBarry Smith           fine_indices[c + 2] = 3 * nidx + 2;
2511e07b27eSBarry Smith           c                   = c + 3;
2521e07b27eSBarry Smith         }
2531e07b27eSBarry Smith       }
2541e07b27eSBarry Smith     }
2551e07b27eSBarry Smith   }
2561e07b27eSBarry Smith 
2571e07b27eSBarry Smith   /* generate scatter */
2589566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), Ml * Nl * Pl * 3, fine_indices, PETSC_USE_POINTER, &is_fine));
2599566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, Ml * Nl * Pl * 3, 0, 1, &is_local));
2601e07b27eSBarry Smith 
2611e07b27eSBarry Smith   /* scatter */
2629566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &perm_coors));
2639566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(perm_coors, PETSC_DECIDE, Ml * Nl * Pl * 3));
2649566063dSJacob Faibussowitsch   PetscCall(VecSetType(perm_coors, VECSEQ));
2659566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(coor_natural, is_fine, perm_coors, is_local, &sctx));
2669566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(sctx, coor_natural, perm_coors, INSERT_VALUES, SCATTER_FORWARD));
2679566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(sctx, coor_natural, perm_coors, INSERT_VALUES, SCATTER_FORWARD));
2681e07b27eSBarry Smith 
2691e07b27eSBarry Smith   /* access */
27057f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
2711e07b27eSBarry Smith     Vec                _coors;
2721e07b27eSBarry Smith     const PetscScalar *LA_perm;
2731e07b27eSBarry Smith     PetscScalar       *LA_coors;
2741e07b27eSBarry Smith 
2759566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(subdm, &_coors));
2769566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(perm_coors, &LA_perm));
2779566063dSJacob Faibussowitsch     PetscCall(VecGetArray(_coors, &LA_coors));
2789371c9d4SSatish Balay     for (i = 0; i < Ml * Nl * Pl * 3; i++) { LA_coors[i] = LA_perm[i]; }
2799566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(_coors, &LA_coors));
2809566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(perm_coors, &LA_perm));
2811e07b27eSBarry Smith   }
2821e07b27eSBarry Smith 
2831e07b27eSBarry Smith   /* update local coords */
28457f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
2851e07b27eSBarry Smith     DM  _dmc;
2861e07b27eSBarry Smith     Vec _coors, _coors_local;
2871e07b27eSBarry Smith 
2889566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(subdm, &_dmc));
2899566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinates(subdm, &_coors));
2909566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(subdm, &_coors_local));
2919566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalBegin(_dmc, _coors, INSERT_VALUES, _coors_local));
2929566063dSJacob Faibussowitsch     PetscCall(DMGlobalToLocalEnd(_dmc, _coors, INSERT_VALUES, _coors_local));
2931e07b27eSBarry Smith   }
2941e07b27eSBarry Smith 
2959566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&sctx));
2969566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_fine));
2979566063dSJacob Faibussowitsch   PetscCall(PetscFree(fine_indices));
2989566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is_local));
2999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&perm_coors));
3009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coor_natural));
3011e07b27eSBarry Smith   PetscFunctionReturn(0);
3021e07b27eSBarry Smith }
3031e07b27eSBarry Smith 
3049371c9d4SSatish Balay static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors(PC pc, PC_Telescope sred, PC_Telescope_DMDACtx *ctx) {
3051e07b27eSBarry Smith   PetscInt     dim;
3061e07b27eSBarry Smith   DM           dm, subdm;
3071e07b27eSBarry Smith   PetscSubcomm psubcomm;
3081e07b27eSBarry Smith   MPI_Comm     comm;
3091e07b27eSBarry Smith   Vec          coor;
3101e07b27eSBarry Smith 
3111e07b27eSBarry Smith   PetscFunctionBegin;
3129566063dSJacob Faibussowitsch   PetscCall(PCGetDM(pc, &dm));
3139566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(dm, &coor));
3141e07b27eSBarry Smith   if (!coor) PetscFunctionReturn(0);
3151e07b27eSBarry Smith   psubcomm = sred->psubcomm;
3161e07b27eSBarry Smith   comm     = PetscSubcommParent(psubcomm);
3171e07b27eSBarry Smith   subdm    = ctx->dmrepart;
3181e07b27eSBarry Smith 
3199566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "PCTelescope: setting up the coordinates (DMDA)\n"));
3209566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, &dim, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
3211e07b27eSBarry Smith   switch (dim) {
3229371c9d4SSatish Balay   case 1: SETERRQ(comm, PETSC_ERR_SUP, "Telescope: DMDA (1D) repartitioning not provided");
3239371c9d4SSatish Balay   case 2: PetscCall(PCTelescopeSetUp_dmda_repart_coors2d(sred, dm, subdm)); break;
3249371c9d4SSatish Balay   case 3: PetscCall(PCTelescopeSetUp_dmda_repart_coors3d(sred, dm, subdm)); break;
3251e07b27eSBarry Smith   }
3261e07b27eSBarry Smith   PetscFunctionReturn(0);
3271e07b27eSBarry Smith }
3281e07b27eSBarry Smith 
3291e07b27eSBarry Smith /* setup repartitioned dm */
3309371c9d4SSatish Balay PetscErrorCode PCTelescopeSetUp_dmda_repart(PC pc, PC_Telescope sred, PC_Telescope_DMDACtx *ctx) {
3311e07b27eSBarry Smith   DM                    dm;
3321e07b27eSBarry Smith   PetscInt              dim, nx, ny, nz, ndof, nsw, sum, k;
3331e07b27eSBarry Smith   DMBoundaryType        bx, by, bz;
3341e07b27eSBarry Smith   DMDAStencilType       stencil;
3351e07b27eSBarry Smith   const PetscInt       *_range_i_re;
3361e07b27eSBarry Smith   const PetscInt       *_range_j_re;
3371e07b27eSBarry Smith   const PetscInt       *_range_k_re;
3381e07b27eSBarry Smith   DMDAInterpolationType itype;
3391e07b27eSBarry Smith   PetscInt              refine_x, refine_y, refine_z;
3401e07b27eSBarry Smith   MPI_Comm              comm, subcomm;
3411e07b27eSBarry Smith   const char           *prefix;
3421e07b27eSBarry Smith 
3431e07b27eSBarry Smith   PetscFunctionBegin;
3441e07b27eSBarry Smith   comm    = PetscSubcommParent(sred->psubcomm);
3451e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
3469566063dSJacob Faibussowitsch   PetscCall(PCGetDM(pc, &dm));
3471e07b27eSBarry Smith 
3489566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, &dim, &nx, &ny, &nz, NULL, NULL, NULL, &ndof, &nsw, &bx, &by, &bz, &stencil));
3499566063dSJacob Faibussowitsch   PetscCall(DMDAGetInterpolationType(dm, &itype));
3509566063dSJacob Faibussowitsch   PetscCall(DMDAGetRefinementFactor(dm, &refine_x, &refine_y, &refine_z));
3511e07b27eSBarry Smith 
3521e07b27eSBarry Smith   ctx->dmrepart = NULL;
3531e07b27eSBarry Smith   _range_i_re = _range_j_re = _range_k_re = NULL;
3541e07b27eSBarry Smith   /* Create DMDA on the child communicator */
35557f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
3561e07b27eSBarry Smith     switch (dim) {
3571e07b27eSBarry Smith     case 1:
3589566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "PCTelescope: setting up the DMDA on comm subset (1D)\n"));
3599566063dSJacob Faibussowitsch       /* PetscCall(DMDACreate1d(subcomm,bx,nx,ndof,nsw,NULL,&ctx->dmrepart)); */
3601e07b27eSBarry Smith       ny = nz = 1;
3611e07b27eSBarry Smith       by = bz = DM_BOUNDARY_NONE;
3621e07b27eSBarry Smith       break;
3631e07b27eSBarry Smith     case 2:
3649566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "PCTelescope: setting up the DMDA on comm subset (2D)\n"));
3659566063dSJacob Faibussowitsch       /* PetscCall(DMDACreate2d(subcomm,bx,by,stencil,nx,ny, PETSC_DECIDE,PETSC_DECIDE,
3669566063dSJacob Faibussowitsch          ndof,nsw, NULL,NULL,&ctx->dmrepart)); */
3671e07b27eSBarry Smith       nz = 1;
3681e07b27eSBarry Smith       bz = DM_BOUNDARY_NONE;
3691e07b27eSBarry Smith       break;
3701e07b27eSBarry Smith     case 3:
3719566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "PCTelescope: setting up the DMDA on comm subset (3D)\n"));
3729566063dSJacob Faibussowitsch       /* PetscCall(DMDACreate3d(subcomm,bx,by,bz,stencil,nx,ny,nz,
3739566063dSJacob Faibussowitsch          PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE, ndof,nsw, NULL,NULL,NULL,&ctx->dmrepart)); */
3741e07b27eSBarry Smith       break;
3751e07b27eSBarry Smith     }
3761e07b27eSBarry Smith     /*
3771e07b27eSBarry Smith      The API DMDACreate1d(), DMDACreate2d(), DMDACreate3d() does not allow us to set/append
3781e07b27eSBarry Smith      a unique option prefix for the DM, thus I prefer to expose the contents of these API's here.
3791e07b27eSBarry Smith      This allows users to control the partitioning of the subDM.
3801e07b27eSBarry Smith     */
3819566063dSJacob Faibussowitsch     PetscCall(DMDACreate(subcomm, &ctx->dmrepart));
3821e07b27eSBarry Smith     /* Set unique option prefix name */
3839566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(sred->ksp, &prefix));
3849566063dSJacob Faibussowitsch     PetscCall(DMSetOptionsPrefix(ctx->dmrepart, prefix));
3859566063dSJacob Faibussowitsch     PetscCall(DMAppendOptionsPrefix(ctx->dmrepart, "repart_"));
3861e07b27eSBarry Smith     /* standard setup from DMDACreate{1,2,3}d() */
3879566063dSJacob Faibussowitsch     PetscCall(DMSetDimension(ctx->dmrepart, dim));
3889566063dSJacob Faibussowitsch     PetscCall(DMDASetSizes(ctx->dmrepart, nx, ny, nz));
3899566063dSJacob Faibussowitsch     PetscCall(DMDASetNumProcs(ctx->dmrepart, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE));
3909566063dSJacob Faibussowitsch     PetscCall(DMDASetBoundaryType(ctx->dmrepart, bx, by, bz));
3919566063dSJacob Faibussowitsch     PetscCall(DMDASetDof(ctx->dmrepart, ndof));
3929566063dSJacob Faibussowitsch     PetscCall(DMDASetStencilType(ctx->dmrepart, stencil));
3939566063dSJacob Faibussowitsch     PetscCall(DMDASetStencilWidth(ctx->dmrepart, nsw));
3949566063dSJacob Faibussowitsch     PetscCall(DMDASetOwnershipRanges(ctx->dmrepart, NULL, NULL, NULL));
3959566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(ctx->dmrepart));
3969566063dSJacob Faibussowitsch     PetscCall(DMSetUp(ctx->dmrepart));
3971e07b27eSBarry Smith     /* Set refinement factors and interpolation type from the partent */
3989566063dSJacob Faibussowitsch     PetscCall(DMDASetRefinementFactor(ctx->dmrepart, refine_x, refine_y, refine_z));
3999566063dSJacob Faibussowitsch     PetscCall(DMDASetInterpolationType(ctx->dmrepart, itype));
4001e07b27eSBarry Smith 
4019566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(ctx->dmrepart, NULL, NULL, NULL, NULL, &ctx->Mp_re, &ctx->Np_re, &ctx->Pp_re, NULL, NULL, NULL, NULL, NULL, NULL));
4029566063dSJacob Faibussowitsch     PetscCall(DMDAGetOwnershipRanges(ctx->dmrepart, &_range_i_re, &_range_j_re, &_range_k_re));
4035e897e82SDave May 
4045e897e82SDave May     ctx->dmrepart->ops->creatematrix              = dm->ops->creatematrix;
4055e897e82SDave May     ctx->dmrepart->ops->createdomaindecomposition = dm->ops->createdomaindecomposition;
4061e07b27eSBarry Smith   }
4071e07b27eSBarry Smith 
4081e07b27eSBarry Smith   /* generate ranges for repartitioned dm */
4091e07b27eSBarry Smith   /* note - assume rank 0 always participates */
410071fcb05SBarry Smith   /* TODO: use a single MPI call */
4119566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&ctx->Mp_re, 1, MPIU_INT, 0, comm));
4129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&ctx->Np_re, 1, MPIU_INT, 0, comm));
4139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&ctx->Pp_re, 1, MPIU_INT, 0, comm));
4141e07b27eSBarry Smith 
4159566063dSJacob Faibussowitsch   PetscCall(PetscCalloc3(ctx->Mp_re, &ctx->range_i_re, ctx->Np_re, &ctx->range_j_re, ctx->Pp_re, &ctx->range_k_re));
4161e07b27eSBarry Smith 
4179566063dSJacob Faibussowitsch   if (_range_i_re) PetscCall(PetscArraycpy(ctx->range_i_re, _range_i_re, ctx->Mp_re));
4189566063dSJacob Faibussowitsch   if (_range_j_re) PetscCall(PetscArraycpy(ctx->range_j_re, _range_j_re, ctx->Np_re));
4199566063dSJacob Faibussowitsch   if (_range_k_re) PetscCall(PetscArraycpy(ctx->range_k_re, _range_k_re, ctx->Pp_re));
4201e07b27eSBarry Smith 
421071fcb05SBarry Smith   /* TODO: use a single MPI call */
4229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(ctx->range_i_re, ctx->Mp_re, MPIU_INT, 0, comm));
4239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(ctx->range_j_re, ctx->Np_re, MPIU_INT, 0, comm));
4249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(ctx->range_k_re, ctx->Pp_re, MPIU_INT, 0, comm));
4251e07b27eSBarry Smith 
4269566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(ctx->Mp_re, &ctx->start_i_re, ctx->Np_re, &ctx->start_j_re, ctx->Pp_re, &ctx->start_k_re));
4271e07b27eSBarry Smith 
4281e07b27eSBarry Smith   sum = 0;
4291e07b27eSBarry Smith   for (k = 0; k < ctx->Mp_re; k++) {
4301e07b27eSBarry Smith     ctx->start_i_re[k] = sum;
4311e07b27eSBarry Smith     sum += ctx->range_i_re[k];
4321e07b27eSBarry Smith   }
4331e07b27eSBarry Smith 
4341e07b27eSBarry Smith   sum = 0;
4351e07b27eSBarry Smith   for (k = 0; k < ctx->Np_re; k++) {
4361e07b27eSBarry Smith     ctx->start_j_re[k] = sum;
4371e07b27eSBarry Smith     sum += ctx->range_j_re[k];
4381e07b27eSBarry Smith   }
4391e07b27eSBarry Smith 
4401e07b27eSBarry Smith   sum = 0;
4411e07b27eSBarry Smith   for (k = 0; k < ctx->Pp_re; k++) {
4421e07b27eSBarry Smith     ctx->start_k_re[k] = sum;
4431e07b27eSBarry Smith     sum += ctx->range_k_re[k];
4441e07b27eSBarry Smith   }
4451e07b27eSBarry Smith 
446ba1c3560SDave May   /* attach repartitioned dm to child ksp */
447ba1c3560SDave May   {
448ba1c3560SDave May     PetscErrorCode (*dmksp_func)(KSP, Mat, Mat, void *);
449ba1c3560SDave May     void *dmksp_ctx;
450ba1c3560SDave May 
4519566063dSJacob Faibussowitsch     PetscCall(DMKSPGetComputeOperators(dm, &dmksp_func, &dmksp_ctx));
452ba1c3560SDave May 
4531e07b27eSBarry Smith     /* attach dm to ksp on sub communicator */
45457f12427SDave May     if (PCTelescope_isActiveRank(sred)) {
4559566063dSJacob Faibussowitsch       PetscCall(KSPSetDM(sred->ksp, ctx->dmrepart));
456ba1c3560SDave May 
457c5db1f53SDave May       if (!dmksp_func || sred->ignore_kspcomputeoperators) {
4589566063dSJacob Faibussowitsch         PetscCall(KSPSetDMActive(sred->ksp, PETSC_FALSE));
459ba1c3560SDave May       } else {
460ba1c3560SDave May         /* sub ksp inherits dmksp_func and context provided by user */
4619566063dSJacob Faibussowitsch         PetscCall(KSPSetComputeOperators(sred->ksp, dmksp_func, dmksp_ctx));
4629566063dSJacob Faibussowitsch         PetscCall(KSPSetDMActive(sred->ksp, PETSC_TRUE));
463ba1c3560SDave May       }
464ba1c3560SDave May     }
4651e07b27eSBarry Smith   }
4661e07b27eSBarry Smith   PetscFunctionReturn(0);
4671e07b27eSBarry Smith }
4681e07b27eSBarry Smith 
4699371c9d4SSatish Balay PetscErrorCode PCTelescopeSetUp_dmda_permutation_3d(PC pc, PC_Telescope sred, PC_Telescope_DMDACtx *ctx) {
4701e07b27eSBarry Smith   DM       dm;
4711e07b27eSBarry Smith   MPI_Comm comm;
4721e07b27eSBarry Smith   Mat      Pscalar, P;
4731e07b27eSBarry Smith   PetscInt ndof;
4741e07b27eSBarry Smith   PetscInt i, j, k, location, startI[3], endI[3], lenI[3], nx, ny, nz;
4751e07b27eSBarry Smith   PetscInt sr, er, Mr;
4761e07b27eSBarry Smith   Vec      V;
4771e07b27eSBarry Smith 
4781e07b27eSBarry Smith   PetscFunctionBegin;
4799566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "PCTelescope: setting up the permutation matrix (DMDA-3D)\n"));
4809566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
4811e07b27eSBarry Smith 
4829566063dSJacob Faibussowitsch   PetscCall(PCGetDM(pc, &dm));
4839566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, NULL, &nx, &ny, &nz, NULL, NULL, NULL, &ndof, NULL, NULL, NULL, NULL, NULL));
4841e07b27eSBarry Smith 
4859566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &V));
4869566063dSJacob Faibussowitsch   PetscCall(VecGetSize(V, &Mr));
4879566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(V, &sr, &er));
4889566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &V));
4891e07b27eSBarry Smith   sr = sr / ndof;
4901e07b27eSBarry Smith   er = er / ndof;
4911e07b27eSBarry Smith   Mr = Mr / ndof;
4921e07b27eSBarry Smith 
4939566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Pscalar));
4949566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Pscalar, (er - sr), (er - sr), Mr, Mr));
4959566063dSJacob Faibussowitsch   PetscCall(MatSetType(Pscalar, MATAIJ));
4969566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Pscalar, 1, NULL));
4979566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Pscalar, 1, NULL, 1, NULL));
4981e07b27eSBarry Smith 
4999566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dm, NULL, NULL, NULL, &lenI[0], &lenI[1], &lenI[2]));
5009566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dm, &startI[0], &startI[1], &startI[2], &endI[0], &endI[1], &endI[2]));
5011e07b27eSBarry Smith   endI[0] += startI[0];
5021e07b27eSBarry Smith   endI[1] += startI[1];
5031e07b27eSBarry Smith   endI[2] += startI[2];
5041e07b27eSBarry Smith 
5051e07b27eSBarry Smith   for (k = startI[2]; k < endI[2]; k++) {
5061e07b27eSBarry Smith     for (j = startI[1]; j < endI[1]; j++) {
5071e07b27eSBarry Smith       for (i = startI[0]; i < endI[0]; i++) {
5081e07b27eSBarry Smith         PetscMPIInt rank_ijk_re, rank_reI[3];
5091e07b27eSBarry Smith         PetscInt    s0_re;
510c6a0d831SBarry Smith         PetscInt    ii, jj, kk, local_ijk_re, mapped_ijk;
5111e07b27eSBarry Smith         PetscInt    lenI_re[3];
5121e07b27eSBarry Smith 
5131e07b27eSBarry Smith         location = (i - startI[0]) + (j - startI[1]) * lenI[0] + (k - startI[2]) * lenI[0] * lenI[1];
5149371c9d4SSatish Balay         PetscCall(_DMDADetermineRankFromGlobalIJK(3, i, j, k, ctx->Mp_re, ctx->Np_re, ctx->Pp_re, ctx->start_i_re, ctx->start_j_re, ctx->start_k_re, ctx->range_i_re, ctx->range_j_re, ctx->range_k_re, &rank_reI[0], &rank_reI[1], &rank_reI[2], &rank_ijk_re));
5159566063dSJacob Faibussowitsch         PetscCall(_DMDADetermineGlobalS0(3, rank_ijk_re, ctx->Mp_re, ctx->Np_re, ctx->Pp_re, ctx->range_i_re, ctx->range_j_re, ctx->range_k_re, &s0_re));
5161e07b27eSBarry Smith         ii = i - ctx->start_i_re[rank_reI[0]];
51708401ef6SPierre Jolivet         PetscCheck(ii >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "[dmdarepart-perm3d] index error ii");
5181e07b27eSBarry Smith         jj = j - ctx->start_j_re[rank_reI[1]];
51908401ef6SPierre Jolivet         PetscCheck(jj >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "[dmdarepart-perm3d] index error jj");
5201e07b27eSBarry Smith         kk = k - ctx->start_k_re[rank_reI[2]];
52108401ef6SPierre Jolivet         PetscCheck(kk >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "[dmdarepart-perm3d] index error kk");
5221e07b27eSBarry Smith         lenI_re[0]   = ctx->range_i_re[rank_reI[0]];
5231e07b27eSBarry Smith         lenI_re[1]   = ctx->range_j_re[rank_reI[1]];
5241e07b27eSBarry Smith         lenI_re[2]   = ctx->range_k_re[rank_reI[2]];
5251e07b27eSBarry Smith         local_ijk_re = ii + jj * lenI_re[0] + kk * lenI_re[0] * lenI_re[1];
5261e07b27eSBarry Smith         mapped_ijk   = s0_re + local_ijk_re;
5279566063dSJacob Faibussowitsch         PetscCall(MatSetValue(Pscalar, sr + location, mapped_ijk, 1.0, INSERT_VALUES));
5281e07b27eSBarry Smith       }
5291e07b27eSBarry Smith     }
5301e07b27eSBarry Smith   }
5319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Pscalar, MAT_FINAL_ASSEMBLY));
5329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Pscalar, MAT_FINAL_ASSEMBLY));
5339566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(Pscalar, ndof, &P));
5349566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Pscalar));
5351e07b27eSBarry Smith   ctx->permutation = P;
5361e07b27eSBarry Smith   PetscFunctionReturn(0);
5371e07b27eSBarry Smith }
5381e07b27eSBarry Smith 
5399371c9d4SSatish Balay PetscErrorCode PCTelescopeSetUp_dmda_permutation_2d(PC pc, PC_Telescope sred, PC_Telescope_DMDACtx *ctx) {
5401e07b27eSBarry Smith   DM       dm;
5411e07b27eSBarry Smith   MPI_Comm comm;
5421e07b27eSBarry Smith   Mat      Pscalar, P;
5431e07b27eSBarry Smith   PetscInt ndof;
5441e07b27eSBarry Smith   PetscInt i, j, location, startI[2], endI[2], lenI[2], nx, ny, nz;
5451e07b27eSBarry Smith   PetscInt sr, er, Mr;
5461e07b27eSBarry Smith   Vec      V;
5471e07b27eSBarry Smith 
5481e07b27eSBarry Smith   PetscFunctionBegin;
5499566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "PCTelescope: setting up the permutation matrix (DMDA-2D)\n"));
5509566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
5519566063dSJacob Faibussowitsch   PetscCall(PCGetDM(pc, &dm));
5529566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, NULL, &nx, &ny, &nz, NULL, NULL, NULL, &ndof, NULL, NULL, NULL, NULL, NULL));
5539566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dm, &V));
5549566063dSJacob Faibussowitsch   PetscCall(VecGetSize(V, &Mr));
5559566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(V, &sr, &er));
5569566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dm, &V));
5571e07b27eSBarry Smith   sr = sr / ndof;
5581e07b27eSBarry Smith   er = er / ndof;
5591e07b27eSBarry Smith   Mr = Mr / ndof;
5601e07b27eSBarry Smith 
5619566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, &Pscalar));
5629566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Pscalar, (er - sr), (er - sr), Mr, Mr));
5639566063dSJacob Faibussowitsch   PetscCall(MatSetType(Pscalar, MATAIJ));
5649566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Pscalar, 1, NULL));
5659566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(Pscalar, 1, NULL, 1, NULL));
5661e07b27eSBarry Smith 
5679566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dm, NULL, NULL, NULL, &lenI[0], &lenI[1], NULL));
5689566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dm, &startI[0], &startI[1], NULL, &endI[0], &endI[1], NULL));
5691e07b27eSBarry Smith   endI[0] += startI[0];
5701e07b27eSBarry Smith   endI[1] += startI[1];
5711e07b27eSBarry Smith 
5721e07b27eSBarry Smith   for (j = startI[1]; j < endI[1]; j++) {
5731e07b27eSBarry Smith     for (i = startI[0]; i < endI[0]; i++) {
5741e07b27eSBarry Smith       PetscMPIInt rank_ijk_re, rank_reI[3];
5751e07b27eSBarry Smith       PetscInt    s0_re;
576c6a0d831SBarry Smith       PetscInt    ii, jj, local_ijk_re, mapped_ijk;
5771e07b27eSBarry Smith       PetscInt    lenI_re[3];
5781e07b27eSBarry Smith 
5791e07b27eSBarry Smith       location = (i - startI[0]) + (j - startI[1]) * lenI[0];
5809371c9d4SSatish Balay       PetscCall(_DMDADetermineRankFromGlobalIJK(2, i, j, 0, ctx->Mp_re, ctx->Np_re, ctx->Pp_re, ctx->start_i_re, ctx->start_j_re, ctx->start_k_re, ctx->range_i_re, ctx->range_j_re, ctx->range_k_re, &rank_reI[0], &rank_reI[1], NULL, &rank_ijk_re));
5811e07b27eSBarry Smith 
5829566063dSJacob Faibussowitsch       PetscCall(_DMDADetermineGlobalS0(2, rank_ijk_re, ctx->Mp_re, ctx->Np_re, ctx->Pp_re, ctx->range_i_re, ctx->range_j_re, ctx->range_k_re, &s0_re));
5831e07b27eSBarry Smith 
5841e07b27eSBarry Smith       ii = i - ctx->start_i_re[rank_reI[0]];
58508401ef6SPierre Jolivet       PetscCheck(ii >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "[dmdarepart-perm2d] index error ii");
5861e07b27eSBarry Smith       jj = j - ctx->start_j_re[rank_reI[1]];
58708401ef6SPierre Jolivet       PetscCheck(jj >= 0, PETSC_COMM_SELF, PETSC_ERR_USER, "[dmdarepart-perm2d] index error jj");
5881e07b27eSBarry Smith 
5891e07b27eSBarry Smith       lenI_re[0]   = ctx->range_i_re[rank_reI[0]];
5901e07b27eSBarry Smith       lenI_re[1]   = ctx->range_j_re[rank_reI[1]];
5911e07b27eSBarry Smith       local_ijk_re = ii + jj * lenI_re[0];
5921e07b27eSBarry Smith       mapped_ijk   = s0_re + local_ijk_re;
5939566063dSJacob Faibussowitsch       PetscCall(MatSetValue(Pscalar, sr + location, mapped_ijk, 1.0, INSERT_VALUES));
5941e07b27eSBarry Smith     }
5951e07b27eSBarry Smith   }
5969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Pscalar, MAT_FINAL_ASSEMBLY));
5979566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Pscalar, MAT_FINAL_ASSEMBLY));
5989566063dSJacob Faibussowitsch   PetscCall(MatCreateMAIJ(Pscalar, ndof, &P));
5999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Pscalar));
6001e07b27eSBarry Smith   ctx->permutation = P;
6011e07b27eSBarry Smith   PetscFunctionReturn(0);
6021e07b27eSBarry Smith }
6031e07b27eSBarry Smith 
6049371c9d4SSatish Balay PetscErrorCode PCTelescopeSetUp_dmda_scatters(PC pc, PC_Telescope sred, PC_Telescope_DMDACtx *ctx) {
6051e07b27eSBarry Smith   Vec        xred, yred, xtmp, x, xp;
6061e07b27eSBarry Smith   VecScatter scatter;
6071e07b27eSBarry Smith   IS         isin;
6081e07b27eSBarry Smith   Mat        B;
6091e07b27eSBarry Smith   PetscInt   m, bs, st, ed;
6101e07b27eSBarry Smith   MPI_Comm   comm;
6111e07b27eSBarry Smith 
6121e07b27eSBarry Smith   PetscFunctionBegin;
6139566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
6149566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &B));
6159566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(B, &x, NULL));
6169566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(B, &bs));
6179566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &xp));
6183ac26c5eSBarry Smith   m    = 0;
6191e07b27eSBarry Smith   xred = NULL;
6201e07b27eSBarry Smith   yred = NULL;
62157f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
6229566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(ctx->dmrepart, &xred));
6239566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(xred, &yred));
6249566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(xred, &st, &ed));
6259566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, ed - st, st, 1, &isin));
6269566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(xred, &m));
6271e07b27eSBarry Smith   } else {
6289566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(x, &st, &ed));
6299566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, 0, st, 1, &isin));
6301e07b27eSBarry Smith   }
6319566063dSJacob Faibussowitsch   PetscCall(ISSetBlockSize(isin, bs));
6329566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &xtmp));
6339566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(xtmp, m, PETSC_DECIDE));
6349566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(xtmp, bs));
6359566063dSJacob Faibussowitsch   PetscCall(VecSetType(xtmp, ((PetscObject)x)->type_name));
6369566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, isin, xtmp, NULL, &scatter));
6371e07b27eSBarry Smith   sred->xred    = xred;
6381e07b27eSBarry Smith   sred->yred    = yred;
6391e07b27eSBarry Smith   sred->isin    = isin;
6401e07b27eSBarry Smith   sred->scatter = scatter;
6411e07b27eSBarry Smith   sred->xtmp    = xtmp;
6421e07b27eSBarry Smith 
6431e07b27eSBarry Smith   ctx->xp = xp;
6449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
6451e07b27eSBarry Smith   PetscFunctionReturn(0);
6461e07b27eSBarry Smith }
6471e07b27eSBarry Smith 
6489371c9d4SSatish Balay PetscErrorCode PCTelescopeSetUp_dmda(PC pc, PC_Telescope sred) {
6491e07b27eSBarry Smith   PC_Telescope_DMDACtx *ctx;
6501e07b27eSBarry Smith   PetscInt              dim;
6511e07b27eSBarry Smith   DM                    dm;
6521e07b27eSBarry Smith   MPI_Comm              comm;
6531e07b27eSBarry Smith 
6541e07b27eSBarry Smith   PetscFunctionBegin;
6559566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "PCTelescope: setup (DMDA)\n"));
6569566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ctx));
6571e07b27eSBarry Smith   sred->dm_ctx = (void *)ctx;
6581e07b27eSBarry Smith 
6599566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
6609566063dSJacob Faibussowitsch   PetscCall(PCGetDM(pc, &dm));
6619566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, &dim, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
6621e07b27eSBarry Smith 
6631e07b27eSBarry Smith   PCTelescopeSetUp_dmda_repart(pc, sred, ctx);
6641e07b27eSBarry Smith   PCTelescopeSetUp_dmda_repart_coors(pc, sred, ctx);
6651e07b27eSBarry Smith   switch (dim) {
6669371c9d4SSatish Balay   case 1: SETERRQ(comm, PETSC_ERR_SUP, "Telescope: DMDA (1D) repartitioning not provided");
6679371c9d4SSatish Balay   case 2: PetscCall(PCTelescopeSetUp_dmda_permutation_2d(pc, sred, ctx)); break;
6689371c9d4SSatish Balay   case 3: PetscCall(PCTelescopeSetUp_dmda_permutation_3d(pc, sred, ctx)); break;
6691e07b27eSBarry Smith   }
6709566063dSJacob Faibussowitsch   PetscCall(PCTelescopeSetUp_dmda_scatters(pc, sred, ctx));
6711e07b27eSBarry Smith   PetscFunctionReturn(0);
6721e07b27eSBarry Smith }
6731e07b27eSBarry Smith 
6749371c9d4SSatish Balay PetscErrorCode PCTelescopeMatCreate_dmda_dmactivefalse(PC pc, PC_Telescope sred, MatReuse reuse, Mat *A) {
6751e07b27eSBarry Smith   PC_Telescope_DMDACtx *ctx;
6761e07b27eSBarry Smith   MPI_Comm              comm, subcomm;
6771e07b27eSBarry Smith   Mat                   Bperm, Bred, B, P;
6781e07b27eSBarry Smith   PetscInt              nr, nc;
6791e07b27eSBarry Smith   IS                    isrow, iscol;
6801e07b27eSBarry Smith   Mat                   Blocal, *_Blocal;
6811e07b27eSBarry Smith 
6821e07b27eSBarry Smith   PetscFunctionBegin;
6839566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "PCTelescope: updating the redundant preconditioned operator (DMDA)\n"));
6849566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
6851e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
6861e07b27eSBarry Smith   ctx     = (PC_Telescope_DMDACtx *)sred->dm_ctx;
6871e07b27eSBarry Smith 
6889566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &B));
6899566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, &nr, &nc));
6901e07b27eSBarry Smith 
6911e07b27eSBarry Smith   P = ctx->permutation;
6929566063dSJacob Faibussowitsch   PetscCall(MatPtAP(B, P, MAT_INITIAL_MATRIX, 1.1, &Bperm));
6931e07b27eSBarry Smith 
6941e07b27eSBarry Smith   /* Get submatrices */
6951e07b27eSBarry Smith   isrow = sred->isin;
6969566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(comm, nc, 0, 1, &iscol));
6971e07b27eSBarry Smith 
6989566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(Bperm, 1, &isrow, &iscol, MAT_INITIAL_MATRIX, &_Blocal));
6991e07b27eSBarry Smith   Blocal = *_Blocal;
7001e07b27eSBarry Smith   Bred   = NULL;
70157f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
7021e07b27eSBarry Smith     PetscInt mm;
7031e07b27eSBarry Smith 
7041e07b27eSBarry Smith     if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; }
7059566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Blocal, &mm, NULL));
7069566063dSJacob Faibussowitsch     /* PetscCall(MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,PETSC_DECIDE,reuse,&Bred)); */
7079566063dSJacob Faibussowitsch     PetscCall(MatCreateMPIMatConcatenateSeqMat(subcomm, Blocal, mm, reuse, &Bred));
7081e07b27eSBarry Smith   }
7091e07b27eSBarry Smith   *A = Bred;
7101e07b27eSBarry Smith 
7119566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
7129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Bperm));
7139566063dSJacob Faibussowitsch   PetscCall(MatDestroyMatrices(1, &_Blocal));
7141e07b27eSBarry Smith   PetscFunctionReturn(0);
7151e07b27eSBarry Smith }
7161e07b27eSBarry Smith 
7179371c9d4SSatish Balay PetscErrorCode PCTelescopeMatCreate_dmda(PC pc, PC_Telescope sred, MatReuse reuse, Mat *A) {
718ba1c3560SDave May   DM dm;
719ba1c3560SDave May   PetscErrorCode (*dmksp_func)(KSP, Mat, Mat, void *);
720ba1c3560SDave May   void *dmksp_ctx;
721ba1c3560SDave May 
722ba1c3560SDave May   PetscFunctionBegin;
7239566063dSJacob Faibussowitsch   PetscCall(PCGetDM(pc, &dm));
7249566063dSJacob Faibussowitsch   PetscCall(DMKSPGetComputeOperators(dm, &dmksp_func, &dmksp_ctx));
725dc9ee9fdSDave May   /* We assume that dmksp_func = NULL, is equivalent to dmActive = PETSC_FALSE */
7267c5279cbSDave May   if (dmksp_func && !sred->ignore_kspcomputeoperators) {
727ba1c3560SDave May     DM  dmrepart;
72828323a89SDave May     Mat Ak;
729ba1c3560SDave May 
730ba1c3560SDave May     *A = NULL;
73157f12427SDave May     if (PCTelescope_isActiveRank(sred)) {
7329566063dSJacob Faibussowitsch       PetscCall(KSPGetDM(sred->ksp, &dmrepart));
733ba1c3560SDave May       if (reuse == MAT_INITIAL_MATRIX) {
7349566063dSJacob Faibussowitsch         PetscCall(DMCreateMatrix(dmrepart, &Ak));
735ba1c3560SDave May         *A = Ak;
736ba1c3560SDave May       } else if (reuse == MAT_REUSE_MATRIX) {
737ba1c3560SDave May         Ak = *A;
738ba1c3560SDave May       }
7395c5dbb1cSDave May       /*
7405c5dbb1cSDave May        There is no need to explicitly assemble the operator now,
7415c5dbb1cSDave May        the sub-KSP will call the method provided to KSPSetComputeOperators() during KSPSetUp()
7425c5dbb1cSDave May       */
743ba1c3560SDave May     }
744ba1c3560SDave May   } else {
7459566063dSJacob Faibussowitsch     PetscCall(PCTelescopeMatCreate_dmda_dmactivefalse(pc, sred, reuse, A));
746ba1c3560SDave May   }
747ba1c3560SDave May   PetscFunctionReturn(0);
748ba1c3560SDave May }
749ba1c3560SDave May 
7509371c9d4SSatish Balay PetscErrorCode PCTelescopeSubNullSpaceCreate_dmda_Telescope(PC pc, PC_Telescope sred, MatNullSpace nullspace, MatNullSpace *sub_nullspace) {
7511e07b27eSBarry Smith   PetscBool             has_const;
752a947c41eSDave May   PetscInt              i, k, n = 0;
7531e07b27eSBarry Smith   const Vec            *vecs;
754c41e779fSDave May   Vec                  *sub_vecs = NULL;
7551e07b27eSBarry Smith   MPI_Comm              subcomm;
7561e07b27eSBarry Smith   PC_Telescope_DMDACtx *ctx;
7571e07b27eSBarry Smith 
7581e07b27eSBarry Smith   PetscFunctionBegin;
7591e07b27eSBarry Smith   ctx     = (PC_Telescope_DMDACtx *)sred->dm_ctx;
7601e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
7619566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceGetVecs(nullspace, &has_const, &n, &vecs));
7621e07b27eSBarry Smith 
76357f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
7641e07b27eSBarry Smith     /* create new vectors */
765*48a46eb9SPierre Jolivet     if (n) PetscCall(VecDuplicateVecs(sred->xred, n, &sub_vecs));
7661e07b27eSBarry Smith   }
7671e07b27eSBarry Smith 
7681e07b27eSBarry Smith   /* copy entries */
7691e07b27eSBarry Smith   for (k = 0; k < n; k++) {
7701e07b27eSBarry Smith     const PetscScalar *x_array;
7711e07b27eSBarry Smith     PetscScalar       *LA_sub_vec;
77213c30530SDave May     PetscInt           st, ed;
7731e07b27eSBarry Smith 
7741e07b27eSBarry Smith     /* permute vector into ordering associated with re-partitioned dmda */
7759566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(ctx->permutation, vecs[k], ctx->xp));
7761e07b27eSBarry Smith 
7771e07b27eSBarry Smith     /* pull in vector x->xtmp */
7789566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(sred->scatter, ctx->xp, sred->xtmp, INSERT_VALUES, SCATTER_FORWARD));
7799566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(sred->scatter, ctx->xp, sred->xtmp, INSERT_VALUES, SCATTER_FORWARD));
7801e07b27eSBarry Smith 
781392968a1SPatrick Sanan     /* copy vector entries into xred */
7829566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(sred->xtmp, &x_array));
783ea2b237eSDave May     if (sub_vecs) {
784ea2b237eSDave May       if (sub_vecs[k]) {
7859566063dSJacob Faibussowitsch         PetscCall(VecGetOwnershipRange(sub_vecs[k], &st, &ed));
7869566063dSJacob Faibussowitsch         PetscCall(VecGetArray(sub_vecs[k], &LA_sub_vec));
7879371c9d4SSatish Balay         for (i = 0; i < ed - st; i++) { LA_sub_vec[i] = x_array[i]; }
7889566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(sub_vecs[k], &LA_sub_vec));
7891e07b27eSBarry Smith       }
790ea2b237eSDave May     }
7919566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(sred->xtmp, &x_array));
7921e07b27eSBarry Smith   }
7931e07b27eSBarry Smith 
79457f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
795d8b9d5b7SPatrick Sanan     /* create new (near) nullspace for redundant object */
7969566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceCreate(subcomm, has_const, n, sub_vecs, sub_nullspace));
7979566063dSJacob Faibussowitsch     PetscCall(VecDestroyVecs(n, &sub_vecs));
79828b400f6SJacob Faibussowitsch     PetscCheck(!nullspace->remove, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope");
79928b400f6SJacob Faibussowitsch     PetscCheck(!nullspace->rmctx, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Propagation of custom remove callback context not supported when propagating (near) nullspaces with PCTelescope");
800d8b9d5b7SPatrick Sanan   }
801392968a1SPatrick Sanan   PetscFunctionReturn(0);
802392968a1SPatrick Sanan }
803392968a1SPatrick Sanan 
8049371c9d4SSatish Balay PetscErrorCode PCTelescopeMatNullSpaceCreate_dmda(PC pc, PC_Telescope sred, Mat sub_mat) {
805392968a1SPatrick Sanan   Mat B;
806392968a1SPatrick Sanan 
807392968a1SPatrick Sanan   PetscFunctionBegin;
8089566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &B));
809392968a1SPatrick Sanan   {
810392968a1SPatrick Sanan     MatNullSpace nullspace, sub_nullspace;
8119566063dSJacob Faibussowitsch     PetscCall(MatGetNullSpace(B, &nullspace));
812392968a1SPatrick Sanan     if (nullspace) {
8139566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "PCTelescope: generating nullspace (DMDA)\n"));
8149566063dSJacob Faibussowitsch       PetscCall(PCTelescopeSubNullSpaceCreate_dmda_Telescope(pc, sred, nullspace, &sub_nullspace));
81557f12427SDave May       if (PCTelescope_isActiveRank(sred)) {
8169566063dSJacob Faibussowitsch         PetscCall(MatSetNullSpace(sub_mat, sub_nullspace));
8179566063dSJacob Faibussowitsch         PetscCall(MatNullSpaceDestroy(&sub_nullspace));
8181e07b27eSBarry Smith       }
819392968a1SPatrick Sanan     }
820392968a1SPatrick Sanan   }
821392968a1SPatrick Sanan   {
822392968a1SPatrick Sanan     MatNullSpace nearnullspace, sub_nearnullspace;
8239566063dSJacob Faibussowitsch     PetscCall(MatGetNearNullSpace(B, &nearnullspace));
824392968a1SPatrick Sanan     if (nearnullspace) {
8259566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "PCTelescope: generating near nullspace (DMDA)\n"));
8269566063dSJacob Faibussowitsch       PetscCall(PCTelescopeSubNullSpaceCreate_dmda_Telescope(pc, sred, nearnullspace, &sub_nearnullspace));
82757f12427SDave May       if (PCTelescope_isActiveRank(sred)) {
8289566063dSJacob Faibussowitsch         PetscCall(MatSetNearNullSpace(sub_mat, sub_nearnullspace));
8299566063dSJacob Faibussowitsch         PetscCall(MatNullSpaceDestroy(&sub_nearnullspace));
830392968a1SPatrick Sanan       }
831392968a1SPatrick Sanan     }
832392968a1SPatrick Sanan   }
8331e07b27eSBarry Smith   PetscFunctionReturn(0);
8341e07b27eSBarry Smith }
8351e07b27eSBarry Smith 
8369371c9d4SSatish Balay PetscErrorCode PCApply_Telescope_dmda(PC pc, Vec x, Vec y) {
8371e07b27eSBarry Smith   PC_Telescope          sred = (PC_Telescope)pc->data;
8381e07b27eSBarry Smith   Mat                   perm;
8391e07b27eSBarry Smith   Vec                   xtmp, xp, xred, yred;
84013c30530SDave May   PetscInt              i, st, ed;
8411e07b27eSBarry Smith   VecScatter            scatter;
8421e07b27eSBarry Smith   PetscScalar          *array;
8431e07b27eSBarry Smith   const PetscScalar    *x_array;
8441e07b27eSBarry Smith   PC_Telescope_DMDACtx *ctx;
8451e07b27eSBarry Smith 
8461e07b27eSBarry Smith   ctx     = (PC_Telescope_DMDACtx *)sred->dm_ctx;
8471e07b27eSBarry Smith   xtmp    = sred->xtmp;
8481e07b27eSBarry Smith   scatter = sred->scatter;
8491e07b27eSBarry Smith   xred    = sred->xred;
8501e07b27eSBarry Smith   yred    = sred->yred;
8511e07b27eSBarry Smith   perm    = ctx->permutation;
8521e07b27eSBarry Smith   xp      = ctx->xp;
8531e07b27eSBarry Smith 
8541e07b27eSBarry Smith   PetscFunctionBegin;
8559566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
85614c9fce5SDave May 
8571e07b27eSBarry Smith   /* permute vector into ordering associated with re-partitioned dmda */
8589566063dSJacob Faibussowitsch   PetscCall(MatMultTranspose(perm, x, xp));
8591e07b27eSBarry Smith 
8601e07b27eSBarry Smith   /* pull in vector x->xtmp */
8619566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, xp, xtmp, INSERT_VALUES, SCATTER_FORWARD));
8629566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, xp, xtmp, INSERT_VALUES, SCATTER_FORWARD));
8631e07b27eSBarry Smith 
864a5b23f4aSJose E. Roman   /* copy vector entries into xred */
8659566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xtmp, &x_array));
8661e07b27eSBarry Smith   if (xred) {
8671e07b27eSBarry Smith     PetscScalar *LA_xred;
8689566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(xred, &st, &ed));
8691e07b27eSBarry Smith 
8709566063dSJacob Faibussowitsch     PetscCall(VecGetArray(xred, &LA_xred));
8719371c9d4SSatish Balay     for (i = 0; i < ed - st; i++) { LA_xred[i] = x_array[i]; }
8729566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(xred, &LA_xred));
8731e07b27eSBarry Smith   }
8749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xtmp, &x_array));
8751e07b27eSBarry Smith 
8761e07b27eSBarry Smith   /* solve */
87757f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
8789566063dSJacob Faibussowitsch     PetscCall(KSPSolve(sred->ksp, xred, yred));
8799566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(sred->ksp, pc, yred));
8801e07b27eSBarry Smith   }
8811e07b27eSBarry Smith 
8821e07b27eSBarry Smith   /* return vector */
8839566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xtmp, &array));
8841e07b27eSBarry Smith   if (yred) {
8851e07b27eSBarry Smith     const PetscScalar *LA_yred;
8869566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(yred, &st, &ed));
8879566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(yred, &LA_yred));
8889371c9d4SSatish Balay     for (i = 0; i < ed - st; i++) { array[i] = LA_yred[i]; }
8899566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(yred, &LA_yred));
8901e07b27eSBarry Smith   }
8919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xtmp, &array));
8929566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, xtmp, xp, INSERT_VALUES, SCATTER_REVERSE));
8939566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, xtmp, xp, INSERT_VALUES, SCATTER_REVERSE));
8949566063dSJacob Faibussowitsch   PetscCall(MatMult(perm, xp, y));
8951e07b27eSBarry Smith   PetscFunctionReturn(0);
8961e07b27eSBarry Smith }
8971e07b27eSBarry Smith 
8989371c9d4SSatish Balay PetscErrorCode PCApplyRichardson_Telescope_dmda(PC pc, Vec x, Vec y, Vec w, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt its, PetscBool zeroguess, PetscInt *outits, PCRichardsonConvergedReason *reason) {
899f650675bSDave May   PC_Telescope          sred = (PC_Telescope)pc->data;
900f650675bSDave May   Mat                   perm;
901a1d91a28SDave May   Vec                   xtmp, xp, yred;
902f650675bSDave May   PetscInt              i, st, ed;
903f650675bSDave May   VecScatter            scatter;
904f650675bSDave May   const PetscScalar    *x_array;
905c41e779fSDave May   PetscBool             default_init_guess_value = PETSC_FALSE;
906f650675bSDave May   PC_Telescope_DMDACtx *ctx;
907f650675bSDave May 
90857f12427SDave May   PetscFunctionBegin;
909f650675bSDave May   ctx     = (PC_Telescope_DMDACtx *)sred->dm_ctx;
910f650675bSDave May   xtmp    = sred->xtmp;
911f650675bSDave May   scatter = sred->scatter;
912f650675bSDave May   yred    = sred->yred;
913f650675bSDave May   perm    = ctx->permutation;
914f650675bSDave May   xp      = ctx->xp;
915f650675bSDave May 
91608401ef6SPierre Jolivet   PetscCheck(its <= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "PCApplyRichardson_Telescope_dmda only supports max_it = 1");
917f650675bSDave May   *reason = (PCRichardsonConvergedReason)0;
918f650675bSDave May 
919f650675bSDave May   if (!zeroguess) {
9209566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "PCTelescopeDMDA: Scattering y for non-zero-initial guess\n"));
921f650675bSDave May     /* permute vector into ordering associated with re-partitioned dmda */
9229566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(perm, y, xp));
923f650675bSDave May 
924f650675bSDave May     /* pull in vector x->xtmp */
9259566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scatter, xp, xtmp, INSERT_VALUES, SCATTER_FORWARD));
9269566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scatter, xp, xtmp, INSERT_VALUES, SCATTER_FORWARD));
927f650675bSDave May 
928a5b23f4aSJose E. Roman     /* copy vector entries into xred */
9299566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(xtmp, &x_array));
930f650675bSDave May     if (yred) {
931f650675bSDave May       PetscScalar *LA_yred;
9329566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(yred, &st, &ed));
9339566063dSJacob Faibussowitsch       PetscCall(VecGetArray(yred, &LA_yred));
9349371c9d4SSatish Balay       for (i = 0; i < ed - st; i++) { LA_yred[i] = x_array[i]; }
9359566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(yred, &LA_yred));
936f650675bSDave May     }
9379566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(xtmp, &x_array));
938f650675bSDave May   }
939f650675bSDave May 
94057f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
9419566063dSJacob Faibussowitsch     PetscCall(KSPGetInitialGuessNonzero(sred->ksp, &default_init_guess_value));
9429566063dSJacob Faibussowitsch     if (!zeroguess) PetscCall(KSPSetInitialGuessNonzero(sred->ksp, PETSC_TRUE));
943f650675bSDave May   }
944f650675bSDave May 
9459566063dSJacob Faibussowitsch   PetscCall(PCApply_Telescope_dmda(pc, x, y));
946f650675bSDave May 
947*48a46eb9SPierre Jolivet   if (PCTelescope_isActiveRank(sred)) PetscCall(KSPSetInitialGuessNonzero(sred->ksp, default_init_guess_value));
948f650675bSDave May 
949f650675bSDave May   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
950f650675bSDave May   *outits = 1;
951f650675bSDave May   PetscFunctionReturn(0);
952f650675bSDave May }
953f650675bSDave May 
9549371c9d4SSatish Balay PetscErrorCode PCReset_Telescope_dmda(PC pc) {
9551e07b27eSBarry Smith   PC_Telescope          sred = (PC_Telescope)pc->data;
9561e07b27eSBarry Smith   PC_Telescope_DMDACtx *ctx;
9571e07b27eSBarry Smith 
9581e07b27eSBarry Smith   PetscFunctionBegin;
9591e07b27eSBarry Smith   ctx = (PC_Telescope_DMDACtx *)sred->dm_ctx;
9609566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx->xp));
9619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->permutation));
9629566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&ctx->dmrepart));
9639566063dSJacob Faibussowitsch   PetscCall(PetscFree3(ctx->range_i_re, ctx->range_j_re, ctx->range_k_re));
9649566063dSJacob Faibussowitsch   PetscCall(PetscFree3(ctx->start_i_re, ctx->start_j_re, ctx->start_k_re));
9651e07b27eSBarry Smith   PetscFunctionReturn(0);
9661e07b27eSBarry Smith }
9671e07b27eSBarry Smith 
9689371c9d4SSatish Balay PetscErrorCode DMView_DA_Short_3d(DM dm, PetscViewer v) {
9691e07b27eSBarry Smith   PetscInt    M, N, P, m, n, p, ndof, nsw;
9701e07b27eSBarry Smith   MPI_Comm    comm;
9711e07b27eSBarry Smith   PetscMPIInt size;
9721e07b27eSBarry Smith   const char *prefix;
9731e07b27eSBarry Smith 
9741e07b27eSBarry Smith   PetscFunctionBegin;
9759566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
9769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
9779566063dSJacob Faibussowitsch   PetscCall(DMGetOptionsPrefix(dm, &prefix));
9789566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, NULL, &M, &N, &P, &m, &n, &p, &ndof, &nsw, NULL, NULL, NULL, NULL));
9799566063dSJacob Faibussowitsch   if (prefix) PetscCall(PetscViewerASCIIPrintf(v, "DMDA Object:    (%s)    %d MPI processes\n", prefix, size));
9809566063dSJacob Faibussowitsch   else PetscCall(PetscViewerASCIIPrintf(v, "DMDA Object:    %d MPI processes\n", size));
98163a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(v, "  M %" PetscInt_FMT " N %" PetscInt_FMT " P %" PetscInt_FMT " m %" PetscInt_FMT " n %" PetscInt_FMT " p %" PetscInt_FMT " dof %" PetscInt_FMT " overlap %" PetscInt_FMT "\n", M, N, P, m, n, p, ndof, nsw));
9821e07b27eSBarry Smith   PetscFunctionReturn(0);
9831e07b27eSBarry Smith }
9841e07b27eSBarry Smith 
9859371c9d4SSatish Balay PetscErrorCode DMView_DA_Short_2d(DM dm, PetscViewer v) {
9861e07b27eSBarry Smith   PetscInt    M, N, m, n, ndof, nsw;
9871e07b27eSBarry Smith   MPI_Comm    comm;
9881e07b27eSBarry Smith   PetscMPIInt size;
9891e07b27eSBarry Smith   const char *prefix;
9901e07b27eSBarry Smith 
9911e07b27eSBarry Smith   PetscFunctionBegin;
9929566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
9939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
9949566063dSJacob Faibussowitsch   PetscCall(DMGetOptionsPrefix(dm, &prefix));
9959566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, NULL, &M, &N, NULL, &m, &n, NULL, &ndof, &nsw, NULL, NULL, NULL, NULL));
9969566063dSJacob Faibussowitsch   if (prefix) PetscCall(PetscViewerASCIIPrintf(v, "DMDA Object:    (%s)    %d MPI processes\n", prefix, size));
9979566063dSJacob Faibussowitsch   else PetscCall(PetscViewerASCIIPrintf(v, "DMDA Object:    %d MPI processes\n", size));
99863a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(v, "  M %" PetscInt_FMT " N %" PetscInt_FMT " m %" PetscInt_FMT " n %" PetscInt_FMT " dof %" PetscInt_FMT " overlap %" PetscInt_FMT "\n", M, N, m, n, ndof, nsw));
9991e07b27eSBarry Smith   PetscFunctionReturn(0);
10001e07b27eSBarry Smith }
10011e07b27eSBarry Smith 
10029371c9d4SSatish Balay PetscErrorCode DMView_DA_Short(DM dm, PetscViewer v) {
10031e07b27eSBarry Smith   PetscInt dim;
10041e07b27eSBarry Smith 
10051e07b27eSBarry Smith   PetscFunctionBegin;
10069566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, &dim, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
10071e07b27eSBarry Smith   switch (dim) {
10089371c9d4SSatish Balay   case 2: PetscCall(DMView_DA_Short_2d(dm, v)); break;
10099371c9d4SSatish Balay   case 3: PetscCall(DMView_DA_Short_3d(dm, v)); break;
10101e07b27eSBarry Smith   }
10111e07b27eSBarry Smith   PetscFunctionReturn(0);
10121e07b27eSBarry Smith }
1013