xref: /petsc/src/ksp/pc/impls/telescope/telescope.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
18d9f7141SDave May #include <petsc/private/petscimpl.h>
2120bdd93SDave May #include <petsc/private/matimpl.h>
36fc41876SBarry Smith #include <petsc/private/pcimpl.h>
41e07b27eSBarry Smith #include <petscksp.h> /*I "petscksp.h" I*/
51e07b27eSBarry Smith #include <petscdm.h>  /*I "petscdm.h" I*/
6575a0592SBarry Smith #include "../src/ksp/pc/impls/telescope/telescope.h"
71e07b27eSBarry Smith 
8bf00f589SPatrick Sanan static PetscBool  cited      = PETSC_FALSE;
99371c9d4SSatish Balay static const char citation[] = "@inproceedings{MaySananRuppKnepleySmith2016,\n"
10bf00f589SPatrick Sanan                                "  title     = {Extreme-Scale Multigrid Components within PETSc},\n"
11bf00f589SPatrick Sanan                                "  author    = {Dave A. May and Patrick Sanan and Karl Rupp and Matthew G. Knepley and Barry F. Smith},\n"
12bf00f589SPatrick Sanan                                "  booktitle = {Proceedings of the Platform for Advanced Scientific Computing Conference},\n"
13bf00f589SPatrick Sanan                                "  series    = {PASC '16},\n"
14bf00f589SPatrick Sanan                                "  isbn      = {978-1-4503-4126-4},\n"
15bf00f589SPatrick Sanan                                "  location  = {Lausanne, Switzerland},\n"
16bf00f589SPatrick Sanan                                "  pages     = {5:1--5:12},\n"
17bf00f589SPatrick Sanan                                "  articleno = {5},\n"
18bf00f589SPatrick Sanan                                "  numpages  = {12},\n"
19a8d69d7bSBarry Smith                                "  url       = {https://doi.acm.org/10.1145/2929908.2929913},\n"
20bf00f589SPatrick Sanan                                "  doi       = {10.1145/2929908.2929913},\n"
21bf00f589SPatrick Sanan                                "  acmid     = {2929913},\n"
22bf00f589SPatrick Sanan                                "  publisher = {ACM},\n"
23bf00f589SPatrick Sanan                                "  address   = {New York, NY, USA},\n"
24bf00f589SPatrick Sanan                                "  keywords  = {GPU, HPC, agglomeration, coarse-level solver, multigrid, parallel computing, preconditioning},\n"
25bf00f589SPatrick Sanan                                "  year      = {2016}\n"
26bf00f589SPatrick Sanan                                "}\n";
27bf00f589SPatrick Sanan 
281e07b27eSBarry Smith /*
29c5083d92SDave May  default setup mode
301e07b27eSBarry Smith 
31c5083d92SDave May  [1a] scatter to (FORWARD)
321e07b27eSBarry Smith  x(comm) -> xtmp(comm)
33c5083d92SDave May  [1b] local copy (to) ranks with color = 0
341e07b27eSBarry Smith  xred(subcomm) <- xtmp
351e07b27eSBarry Smith 
36c5083d92SDave May  [2] solve on sub KSP to obtain yred(subcomm)
37c5083d92SDave May 
38c5083d92SDave May  [3a] local copy (from) ranks with color = 0
391e07b27eSBarry Smith  yred(subcomm) --> xtmp
40c5083d92SDave May  [2b] scatter from (REVERSE)
411e07b27eSBarry Smith  xtmp(comm) -> y(comm)
421e07b27eSBarry Smith */
431e07b27eSBarry Smith 
448d9f7141SDave May /*
45d083f849SBarry Smith   Collective[comm_f]
468d9f7141SDave May   Notes
478d9f7141SDave May    * Using comm_f = MPI_COMM_NULL will result in an error
488d9f7141SDave May    * Using comm_c = MPI_COMM_NULL is valid. If all instances of comm_c are NULL the subcomm is not valid.
498d9f7141SDave May    * If any non NULL comm_c communicator cannot map any of its ranks to comm_f, the subcomm is not valid.
508d9f7141SDave May */
519371c9d4SSatish Balay PetscErrorCode PCTelescopeTestValidSubcomm(MPI_Comm comm_f, MPI_Comm comm_c, PetscBool *isvalid) {
5257f12427SDave May   PetscInt     valid = 1;
538d9f7141SDave May   MPI_Group    group_f, group_c;
548d9f7141SDave May   PetscMPIInt  count, k, size_f = 0, size_c = 0, size_c_sum = 0;
555bd1e576SStefano Zampini   PetscMPIInt *ranks_f, *ranks_c;
568d9f7141SDave May 
5757f12427SDave May   PetscFunctionBegin;
5808401ef6SPierre Jolivet   PetscCheck(comm_f != MPI_COMM_NULL, PETSC_COMM_SELF, PETSC_ERR_SUP, "comm_f cannot be MPI_COMM_NULL");
598d9f7141SDave May 
609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_group(comm_f, &group_f));
61*48a46eb9SPierre Jolivet   if (comm_c != MPI_COMM_NULL) PetscCallMPI(MPI_Comm_group(comm_c, &group_c));
628d9f7141SDave May 
639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm_f, &size_f));
64*48a46eb9SPierre Jolivet   if (comm_c != MPI_COMM_NULL) PetscCallMPI(MPI_Comm_size(comm_c, &size_c));
658d9f7141SDave May 
668d9f7141SDave May   /* check not all comm_c's are NULL */
678d9f7141SDave May   size_c_sum = size_c;
689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &size_c_sum, 1, MPI_INT, MPI_SUM, comm_f));
695bd1e576SStefano Zampini   if (size_c_sum == 0) valid = 0;
708d9f7141SDave May 
718d9f7141SDave May   /* check we can map at least 1 rank in comm_c to comm_f */
729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size_f, &ranks_f));
739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size_c, &ranks_c));
745bd1e576SStefano Zampini   for (k = 0; k < size_f; k++) ranks_f[k] = MPI_UNDEFINED;
755bd1e576SStefano Zampini   for (k = 0; k < size_c; k++) ranks_c[k] = k;
768d9f7141SDave May 
7757f12427SDave May   /*
7857f12427SDave May    MPI_Group_translate_ranks() returns a non-zero exit code if any rank cannot be translated.
7957f12427SDave May    I do not want the code to terminate immediately if this occurs, rather I want to throw
8057f12427SDave May    the error later (during PCSetUp_Telescope()) via SETERRQ() with a message indicating
8157f12427SDave May    that comm_c is not a valid sub-communicator.
829566063dSJacob Faibussowitsch    Hence I purposefully do not call PetscCall() after MPI_Group_translate_ranks().
8357f12427SDave May   */
848d9f7141SDave May   count = 0;
858d9f7141SDave May   if (comm_c != MPI_COMM_NULL) {
8666b79024SDave May     (void)MPI_Group_translate_ranks(group_c, size_c, ranks_c, group_f, ranks_f);
878d9f7141SDave May     for (k = 0; k < size_f; k++) {
889371c9d4SSatish Balay       if (ranks_f[k] == MPI_UNDEFINED) { count++; }
898d9f7141SDave May     }
908d9f7141SDave May   }
915bd1e576SStefano Zampini   if (count == size_f) valid = 0;
928d9f7141SDave May 
939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &valid, 1, MPIU_INT, MPI_MIN, comm_f));
945bd1e576SStefano Zampini   if (valid == 1) *isvalid = PETSC_TRUE;
955bd1e576SStefano Zampini   else *isvalid = PETSC_FALSE;
968d9f7141SDave May 
979566063dSJacob Faibussowitsch   PetscCall(PetscFree(ranks_f));
989566063dSJacob Faibussowitsch   PetscCall(PetscFree(ranks_c));
999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Group_free(&group_f));
100*48a46eb9SPierre Jolivet   if (comm_c != MPI_COMM_NULL) PetscCallMPI(MPI_Group_free(&group_c));
1018d9f7141SDave May   PetscFunctionReturn(0);
1028d9f7141SDave May }
1038d9f7141SDave May 
1049371c9d4SSatish Balay DM private_PCTelescopeGetSubDM(PC_Telescope sred) {
105c6a0d831SBarry Smith   DM subdm = NULL;
1061e07b27eSBarry Smith 
1079371c9d4SSatish Balay   if (!PCTelescope_isActiveRank(sred)) {
1089371c9d4SSatish Balay     subdm = NULL;
1099371c9d4SSatish Balay   } else {
1101e07b27eSBarry Smith     switch (sred->sr_type) {
1119371c9d4SSatish Balay     case TELESCOPE_DEFAULT: subdm = NULL; break;
1129371c9d4SSatish Balay     case TELESCOPE_DMDA: subdm = ((PC_Telescope_DMDACtx *)sred->dm_ctx)->dmrepart; break;
1139371c9d4SSatish Balay     case TELESCOPE_DMPLEX: subdm = NULL; break;
1149371c9d4SSatish Balay     case TELESCOPE_COARSEDM:
1159371c9d4SSatish Balay       if (sred->ksp) { KSPGetDM(sred->ksp, &subdm); }
1168d9f7141SDave May       break;
1171e07b27eSBarry Smith     }
1181e07b27eSBarry Smith   }
1191e07b27eSBarry Smith   return (subdm);
1201e07b27eSBarry Smith }
1211e07b27eSBarry Smith 
1229371c9d4SSatish Balay PetscErrorCode PCTelescopeSetUp_default(PC pc, PC_Telescope sred) {
1231e07b27eSBarry Smith   PetscInt   m, M, bs, st, ed;
1241e07b27eSBarry Smith   Vec        x, xred, yred, xtmp;
1251e07b27eSBarry Smith   Mat        B;
1261e07b27eSBarry Smith   MPI_Comm   comm, subcomm;
1271e07b27eSBarry Smith   VecScatter scatter;
1281e07b27eSBarry Smith   IS         isin;
12954cd43dcSJunchao Zhang   VecType    vectype;
1301e07b27eSBarry Smith 
1311e07b27eSBarry Smith   PetscFunctionBegin;
1329566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "PCTelescope: setup (default)\n"));
1331e07b27eSBarry Smith   comm    = PetscSubcommParent(sred->psubcomm);
1341e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
1351e07b27eSBarry Smith 
1369566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &B));
1379566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, &M, NULL));
1389566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(B, &bs));
1399566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(B, &x, NULL));
1409566063dSJacob Faibussowitsch   PetscCall(MatGetVecType(B, &vectype));
1411e07b27eSBarry Smith 
1421e07b27eSBarry Smith   xred = NULL;
1433ac26c5eSBarry Smith   m    = 0;
14457f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
1459566063dSJacob Faibussowitsch     PetscCall(VecCreate(subcomm, &xred));
1469566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(xred, PETSC_DECIDE, M));
1479566063dSJacob Faibussowitsch     PetscCall(VecSetBlockSize(xred, bs));
1489566063dSJacob Faibussowitsch     PetscCall(VecSetType(xred, vectype)); /* Use the preconditioner matrix's vectype by default */
1499566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(xred));
1509566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(xred, &m));
1511e07b27eSBarry Smith   }
1521e07b27eSBarry Smith 
1531e07b27eSBarry Smith   yred = NULL;
154*48a46eb9SPierre Jolivet   if (PCTelescope_isActiveRank(sred)) PetscCall(VecDuplicate(xred, &yred));
1551e07b27eSBarry Smith 
1569566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &xtmp));
1579566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(xtmp, m, PETSC_DECIDE));
1589566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(xtmp, bs));
1599566063dSJacob Faibussowitsch   PetscCall(VecSetType(xtmp, vectype));
1601e07b27eSBarry Smith 
16157f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
1629566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(xred, &st, &ed));
1639566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, (ed - st), st, 1, &isin));
1641e07b27eSBarry Smith   } else {
1659566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(x, &st, &ed));
1669566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, 0, st, 1, &isin));
1671e07b27eSBarry Smith   }
1689566063dSJacob Faibussowitsch   PetscCall(ISSetBlockSize(isin, bs));
1691e07b27eSBarry Smith 
1709566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, isin, xtmp, NULL, &scatter));
1711e07b27eSBarry Smith 
1721e07b27eSBarry Smith   sred->isin    = isin;
1731e07b27eSBarry Smith   sred->scatter = scatter;
1741e07b27eSBarry Smith   sred->xred    = xred;
1751e07b27eSBarry Smith   sred->yred    = yred;
1761e07b27eSBarry Smith   sred->xtmp    = xtmp;
1779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1781e07b27eSBarry Smith   PetscFunctionReturn(0);
1791e07b27eSBarry Smith }
1801e07b27eSBarry Smith 
1819371c9d4SSatish Balay PetscErrorCode PCTelescopeMatCreate_default(PC pc, PC_Telescope sred, MatReuse reuse, Mat *A) {
1821e07b27eSBarry Smith   MPI_Comm comm, subcomm;
1831e07b27eSBarry Smith   Mat      Bred, B;
1845624f943SPierre Jolivet   PetscInt nr, nc, bs;
1851e07b27eSBarry Smith   IS       isrow, iscol;
1861e07b27eSBarry Smith   Mat      Blocal, *_Blocal;
1871e07b27eSBarry Smith 
1881e07b27eSBarry Smith   PetscFunctionBegin;
1899566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "PCTelescope: updating the redundant preconditioned operator (default)\n"));
1909566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
1911e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
1929566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &B));
1939566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, &nr, &nc));
1941e07b27eSBarry Smith   isrow = sred->isin;
1959566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, nc, 0, 1, &iscol));
1969566063dSJacob Faibussowitsch   PetscCall(ISSetIdentity(iscol));
1979566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(B, NULL, &bs));
1989566063dSJacob Faibussowitsch   PetscCall(ISSetBlockSize(iscol, bs));
1999566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_SUBMAT_SINGLEIS, PETSC_TRUE));
2009566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, 1, &isrow, &iscol, MAT_INITIAL_MATRIX, &_Blocal));
2011e07b27eSBarry Smith   Blocal = *_Blocal;
2029566063dSJacob Faibussowitsch   PetscCall(PetscFree(_Blocal));
2031e07b27eSBarry Smith   Bred = NULL;
20457f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
2051e07b27eSBarry Smith     PetscInt mm;
2061e07b27eSBarry Smith 
2071e07b27eSBarry Smith     if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; }
2081e07b27eSBarry Smith 
2099566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Blocal, &mm, NULL));
2109566063dSJacob Faibussowitsch     PetscCall(MatCreateMPIMatConcatenateSeqMat(subcomm, Blocal, mm, reuse, &Bred));
2111e07b27eSBarry Smith   }
2121e07b27eSBarry Smith   *A = Bred;
2139566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
2149566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Blocal));
2151e07b27eSBarry Smith   PetscFunctionReturn(0);
2161e07b27eSBarry Smith }
2171e07b27eSBarry Smith 
2189371c9d4SSatish Balay static PetscErrorCode PCTelescopeSubNullSpaceCreate_Telescope(PC pc, PC_Telescope sred, MatNullSpace nullspace, MatNullSpace *sub_nullspace) {
2191e07b27eSBarry Smith   PetscBool  has_const;
2201e07b27eSBarry Smith   const Vec *vecs;
221c41e779fSDave May   Vec       *sub_vecs = NULL;
222392968a1SPatrick Sanan   PetscInt   i, k, n = 0;
2231e07b27eSBarry Smith   MPI_Comm   subcomm;
2241e07b27eSBarry Smith 
2251e07b27eSBarry Smith   PetscFunctionBegin;
2261e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
2279566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceGetVecs(nullspace, &has_const, &n, &vecs));
2281e07b27eSBarry Smith 
22957f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
230*48a46eb9SPierre Jolivet     if (n) PetscCall(VecDuplicateVecs(sred->xred, n, &sub_vecs));
2311e07b27eSBarry Smith   }
2321e07b27eSBarry Smith 
2331e07b27eSBarry Smith   /* copy entries */
2341e07b27eSBarry Smith   for (k = 0; k < n; k++) {
2351e07b27eSBarry Smith     const PetscScalar *x_array;
2361e07b27eSBarry Smith     PetscScalar       *LA_sub_vec;
23713c30530SDave May     PetscInt           st, ed;
2381e07b27eSBarry Smith 
2391e07b27eSBarry Smith     /* pull in vector x->xtmp */
2409566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(sred->scatter, vecs[k], sred->xtmp, INSERT_VALUES, SCATTER_FORWARD));
2419566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(sred->scatter, vecs[k], sred->xtmp, INSERT_VALUES, SCATTER_FORWARD));
24247856c66SBarry Smith     if (sub_vecs) {
243a04a6428SPatrick Sanan       /* copy vector entries into xred */
2449566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(sred->xtmp, &x_array));
245ea2b237eSDave May       if (sub_vecs[k]) {
2469566063dSJacob Faibussowitsch         PetscCall(VecGetOwnershipRange(sub_vecs[k], &st, &ed));
2479566063dSJacob Faibussowitsch         PetscCall(VecGetArray(sub_vecs[k], &LA_sub_vec));
2489371c9d4SSatish Balay         for (i = 0; i < ed - st; i++) { LA_sub_vec[i] = x_array[i]; }
2499566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(sub_vecs[k], &LA_sub_vec));
2501e07b27eSBarry Smith       }
2519566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(sred->xtmp, &x_array));
2521e07b27eSBarry Smith     }
25347856c66SBarry Smith   }
2541e07b27eSBarry Smith 
25557f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
256d8b9d5b7SPatrick Sanan     /* create new (near) nullspace for redundant object */
2579566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceCreate(subcomm, has_const, n, sub_vecs, sub_nullspace));
2589566063dSJacob Faibussowitsch     PetscCall(VecDestroyVecs(n, &sub_vecs));
25928b400f6SJacob Faibussowitsch     PetscCheck(!nullspace->remove, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope");
26028b400f6SJacob Faibussowitsch     PetscCheck(!nullspace->rmctx, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Propagation of custom remove callback context not supported when propagating (near) nullspaces with PCTelescope");
261d8b9d5b7SPatrick Sanan   }
262392968a1SPatrick Sanan   PetscFunctionReturn(0);
263392968a1SPatrick Sanan }
264392968a1SPatrick Sanan 
2659371c9d4SSatish Balay static PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc, PC_Telescope sred, Mat sub_mat) {
266392968a1SPatrick Sanan   Mat B;
267392968a1SPatrick Sanan 
268392968a1SPatrick Sanan   PetscFunctionBegin;
2699566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &B));
270392968a1SPatrick Sanan   /* Propagate the nullspace if it exists */
271392968a1SPatrick Sanan   {
272392968a1SPatrick Sanan     MatNullSpace nullspace, sub_nullspace;
2739566063dSJacob Faibussowitsch     PetscCall(MatGetNullSpace(B, &nullspace));
274392968a1SPatrick Sanan     if (nullspace) {
2759566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "PCTelescope: generating nullspace (default)\n"));
2769566063dSJacob Faibussowitsch       PetscCall(PCTelescopeSubNullSpaceCreate_Telescope(pc, sred, nullspace, &sub_nullspace));
27757f12427SDave May       if (PCTelescope_isActiveRank(sred)) {
2789566063dSJacob Faibussowitsch         PetscCall(MatSetNullSpace(sub_mat, sub_nullspace));
2799566063dSJacob Faibussowitsch         PetscCall(MatNullSpaceDestroy(&sub_nullspace));
2801e07b27eSBarry Smith       }
281392968a1SPatrick Sanan     }
282392968a1SPatrick Sanan   }
283392968a1SPatrick Sanan   /* Propagate the near nullspace if it exists */
284392968a1SPatrick Sanan   {
285392968a1SPatrick Sanan     MatNullSpace nearnullspace, sub_nearnullspace;
2869566063dSJacob Faibussowitsch     PetscCall(MatGetNearNullSpace(B, &nearnullspace));
287392968a1SPatrick Sanan     if (nearnullspace) {
2889566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "PCTelescope: generating near nullspace (default)\n"));
2899566063dSJacob Faibussowitsch       PetscCall(PCTelescopeSubNullSpaceCreate_Telescope(pc, sred, nearnullspace, &sub_nearnullspace));
29057f12427SDave May       if (PCTelescope_isActiveRank(sred)) {
2919566063dSJacob Faibussowitsch         PetscCall(MatSetNearNullSpace(sub_mat, sub_nearnullspace));
2929566063dSJacob Faibussowitsch         PetscCall(MatNullSpaceDestroy(&sub_nearnullspace));
293392968a1SPatrick Sanan       }
294392968a1SPatrick Sanan     }
295392968a1SPatrick Sanan   }
2961e07b27eSBarry Smith   PetscFunctionReturn(0);
2971e07b27eSBarry Smith }
2981e07b27eSBarry Smith 
2999371c9d4SSatish Balay static PetscErrorCode PCView_Telescope(PC pc, PetscViewer viewer) {
3001e07b27eSBarry Smith   PC_Telescope sred = (PC_Telescope)pc->data;
3011e07b27eSBarry Smith   PetscBool    iascii, isstring;
3021e07b27eSBarry Smith   PetscViewer  subviewer;
3031e07b27eSBarry Smith 
3041e07b27eSBarry Smith   PetscFunctionBegin;
3059566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
3069566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
3071e07b27eSBarry Smith   if (iascii) {
3088d9f7141SDave May     {
3091e07b27eSBarry Smith       MPI_Comm    comm, subcomm;
3101e07b27eSBarry Smith       PetscMPIInt comm_size, subcomm_size;
3118d9f7141SDave May       DM          dm = NULL, subdm = NULL;
3121e07b27eSBarry Smith 
3139566063dSJacob Faibussowitsch       PetscCall(PCGetDM(pc, &dm));
3141e07b27eSBarry Smith       subdm = private_PCTelescopeGetSubDM(sred);
3158d9f7141SDave May 
3168d9f7141SDave May       if (sred->psubcomm) {
3171e07b27eSBarry Smith         comm    = PetscSubcommParent(sred->psubcomm);
3181e07b27eSBarry Smith         subcomm = PetscSubcommChild(sred->psubcomm);
3199566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(comm, &comm_size));
3209566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(subcomm, &subcomm_size));
3211e07b27eSBarry Smith 
3229566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
32363a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "petsc subcomm: parent comm size reduction factor = %" PetscInt_FMT "\n", sred->redfactor));
3249566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "petsc subcomm: parent_size = %d , subcomm_size = %d\n", (int)comm_size, (int)subcomm_size));
32548a10b22SPatrick Sanan         switch (sred->subcommtype) {
3269371c9d4SSatish Balay         case PETSC_SUBCOMM_INTERLACED: PetscCall(PetscViewerASCIIPrintf(viewer, "petsc subcomm: type = %s\n", PetscSubcommTypes[sred->subcommtype])); break;
3279371c9d4SSatish Balay         case PETSC_SUBCOMM_CONTIGUOUS: PetscCall(PetscViewerASCIIPrintf(viewer, "petsc subcomm type = %s\n", PetscSubcommTypes[sred->subcommtype])); break;
3289371c9d4SSatish Balay         default: SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "General subcomm type not supported by PCTelescope");
32948a10b22SPatrick Sanan         }
3309566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
3318d9f7141SDave May       } else {
3329566063dSJacob Faibussowitsch         PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
3338d9f7141SDave May         subcomm = sred->subcomm;
3349371c9d4SSatish Balay         if (!PCTelescope_isActiveRank(sred)) { subcomm = PETSC_COMM_SELF; }
3358d9f7141SDave May 
3369566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
3379566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "subcomm: using user provided sub-communicator\n"));
3389566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
3398d9f7141SDave May       }
3408d9f7141SDave May 
3419566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, subcomm, &subviewer));
34257f12427SDave May       if (PCTelescope_isActiveRank(sred)) {
3439566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(subviewer));
3441e07b27eSBarry Smith 
345*48a46eb9SPierre Jolivet         if (dm && sred->ignore_dm) PetscCall(PetscViewerASCIIPrintf(subviewer, "ignoring DM\n"));
346*48a46eb9SPierre Jolivet         if (sred->ignore_kspcomputeoperators) PetscCall(PetscViewerASCIIPrintf(subviewer, "ignoring KSPComputeOperators\n"));
3471e07b27eSBarry Smith         switch (sred->sr_type) {
3489371c9d4SSatish Balay         case TELESCOPE_DEFAULT: PetscCall(PetscViewerASCIIPrintf(subviewer, "setup type: default\n")); break;
3491e07b27eSBarry Smith         case TELESCOPE_DMDA:
3509566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "setup type: DMDA auto-repartitioning\n"));
3519566063dSJacob Faibussowitsch           PetscCall(DMView_DA_Short(subdm, subviewer));
3521e07b27eSBarry Smith           break;
3539371c9d4SSatish Balay         case TELESCOPE_DMPLEX: PetscCall(PetscViewerASCIIPrintf(subviewer, "setup type: DMPLEX auto-repartitioning\n")); break;
3549371c9d4SSatish Balay         case TELESCOPE_COARSEDM: PetscCall(PetscViewerASCIIPrintf(subviewer, "setup type: coarse DM\n")); break;
3558d9f7141SDave May         }
3568d9f7141SDave May 
3578d9f7141SDave May         if (dm) {
3588d9f7141SDave May           PetscObject obj = (PetscObject)dm;
3599566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "Parent DM object:"));
3608d9f7141SDave May           PetscViewerASCIIUseTabs(subviewer, PETSC_FALSE);
3618d9f7141SDave May           if (obj->type_name) { PetscViewerASCIIPrintf(subviewer, " type = %s;", obj->type_name); }
3628d9f7141SDave May           if (obj->name) { PetscViewerASCIIPrintf(subviewer, " name = %s;", obj->name); }
3638d9f7141SDave May           if (obj->prefix) { PetscViewerASCIIPrintf(subviewer, " prefix = %s", obj->prefix); }
3649566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "\n"));
3658d9f7141SDave May           PetscViewerASCIIUseTabs(subviewer, PETSC_TRUE);
3668d9f7141SDave May         } else {
3679566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "Parent DM object: NULL\n"));
3688d9f7141SDave May         }
3698d9f7141SDave May         if (subdm) {
3708d9f7141SDave May           PetscObject obj = (PetscObject)subdm;
3719566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "Sub DM object:"));
3728d9f7141SDave May           PetscViewerASCIIUseTabs(subviewer, PETSC_FALSE);
3738d9f7141SDave May           if (obj->type_name) { PetscViewerASCIIPrintf(subviewer, " type = %s;", obj->type_name); }
3748d9f7141SDave May           if (obj->name) { PetscViewerASCIIPrintf(subviewer, " name = %s;", obj->name); }
3758d9f7141SDave May           if (obj->prefix) { PetscViewerASCIIPrintf(subviewer, " prefix = %s", obj->prefix); }
3769566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "\n"));
3778d9f7141SDave May           PetscViewerASCIIUseTabs(subviewer, PETSC_TRUE);
3788d9f7141SDave May         } else {
3799566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "Sub DM object: NULL\n"));
3801e07b27eSBarry Smith         }
3811e07b27eSBarry Smith 
3829566063dSJacob Faibussowitsch         PetscCall(KSPView(sred->ksp, subviewer));
3839566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(subviewer));
3841e07b27eSBarry Smith       }
3859566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, subcomm, &subviewer));
3861e07b27eSBarry Smith     }
3871e07b27eSBarry Smith   }
3881e07b27eSBarry Smith   PetscFunctionReturn(0);
3891e07b27eSBarry Smith }
3901e07b27eSBarry Smith 
3919371c9d4SSatish Balay static PetscErrorCode PCSetUp_Telescope(PC pc) {
3921e07b27eSBarry Smith   PC_Telescope    sred = (PC_Telescope)pc->data;
393bd49479cSSatish Balay   MPI_Comm        comm, subcomm = 0;
3941e07b27eSBarry Smith   PCTelescopeType sr_type;
3951e07b27eSBarry Smith 
3961e07b27eSBarry Smith   PetscFunctionBegin;
3979566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
3981e07b27eSBarry Smith 
3991e07b27eSBarry Smith   /* Determine type of setup/update */
4001e07b27eSBarry Smith   if (!pc->setupcalled) {
4011e07b27eSBarry Smith     PetscBool has_dm, same;
4021e07b27eSBarry Smith     DM        dm;
4031e07b27eSBarry Smith 
4041e07b27eSBarry Smith     sr_type = TELESCOPE_DEFAULT;
4051e07b27eSBarry Smith     has_dm  = PETSC_FALSE;
4069566063dSJacob Faibussowitsch     PetscCall(PCGetDM(pc, &dm));
4071e07b27eSBarry Smith     if (dm) { has_dm = PETSC_TRUE; }
4081e07b27eSBarry Smith     if (has_dm) {
4091e07b27eSBarry Smith       /* check for dmda */
4109566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &same));
4111e07b27eSBarry Smith       if (same) {
4129566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "PCTelescope: found DMDA\n"));
4131e07b27eSBarry Smith         sr_type = TELESCOPE_DMDA;
4141e07b27eSBarry Smith       }
4151e07b27eSBarry Smith       /* check for dmplex */
4169566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &same));
4171e07b27eSBarry Smith       if (same) {
4189566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "PCTelescope: found DMPLEX\n"));
4191e07b27eSBarry Smith         sr_type = TELESCOPE_DMPLEX;
4201e07b27eSBarry Smith       }
4218d9f7141SDave May 
4228d9f7141SDave May       if (sred->use_coarse_dm) {
4239566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "PCTelescope: using coarse DM\n"));
4248d9f7141SDave May         sr_type = TELESCOPE_COARSEDM;
4251e07b27eSBarry Smith       }
4261e07b27eSBarry Smith 
4271e07b27eSBarry Smith       if (sred->ignore_dm) {
4289566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "PCTelescope: ignoring DM\n"));
4291e07b27eSBarry Smith         sr_type = TELESCOPE_DEFAULT;
4301e07b27eSBarry Smith       }
4318d9f7141SDave May     }
4321e07b27eSBarry Smith     sred->sr_type = sr_type;
4331e07b27eSBarry Smith   } else {
4341e07b27eSBarry Smith     sr_type = sred->sr_type;
4351e07b27eSBarry Smith   }
4361e07b27eSBarry Smith 
437d8b9d5b7SPatrick Sanan   /* set function pointers for repartition setup, matrix creation/update, matrix (near) nullspace, and reset functionality */
4381e07b27eSBarry Smith   switch (sr_type) {
4391e07b27eSBarry Smith   case TELESCOPE_DEFAULT:
4401e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
4411e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
4421e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
4431e07b27eSBarry Smith     sred->pctelescope_reset_type              = NULL;
4441e07b27eSBarry Smith     break;
4451e07b27eSBarry Smith   case TELESCOPE_DMDA:
4461e07b27eSBarry Smith     pc->ops->apply                            = PCApply_Telescope_dmda;
447f650675bSDave May     pc->ops->applyrichardson                  = PCApplyRichardson_Telescope_dmda;
4481e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_dmda;
4491e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_dmda;
4501e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda;
4511e07b27eSBarry Smith     sred->pctelescope_reset_type              = PCReset_Telescope_dmda;
4521e07b27eSBarry Smith     break;
4539371c9d4SSatish Balay   case TELESCOPE_DMPLEX: SETERRQ(comm, PETSC_ERR_SUP, "Support for DMPLEX is currently not available");
4548d9f7141SDave May   case TELESCOPE_COARSEDM:
4558d9f7141SDave May     pc->ops->apply                            = PCApply_Telescope_CoarseDM;
4568d9f7141SDave May     pc->ops->applyrichardson                  = PCApplyRichardson_Telescope_CoarseDM;
4578d9f7141SDave May     sred->pctelescope_setup_type              = PCTelescopeSetUp_CoarseDM;
4588d9f7141SDave May     sred->pctelescope_matcreate_type          = NULL;
4598d9f7141SDave May     sred->pctelescope_matnullspacecreate_type = NULL; /* PCTelescopeMatNullSpaceCreate_CoarseDM; */
4608d9f7141SDave May     sred->pctelescope_reset_type              = PCReset_Telescope_CoarseDM;
4611e07b27eSBarry Smith     break;
4629371c9d4SSatish Balay   default: SETERRQ(comm, PETSC_ERR_SUP, "Support only provided for: repartitioning an operator; repartitioning a DMDA; or using a coarse DM");
4638d9f7141SDave May   }
4648d9f7141SDave May 
4658d9f7141SDave May   /* subcomm definition */
4668d9f7141SDave May   if (!pc->setupcalled) {
4678d9f7141SDave May     if ((sr_type == TELESCOPE_DEFAULT) || (sr_type == TELESCOPE_DMDA)) {
4688d9f7141SDave May       if (!sred->psubcomm) {
4699566063dSJacob Faibussowitsch         PetscCall(PetscSubcommCreate(comm, &sred->psubcomm));
4709566063dSJacob Faibussowitsch         PetscCall(PetscSubcommSetNumber(sred->psubcomm, sred->redfactor));
4719566063dSJacob Faibussowitsch         PetscCall(PetscSubcommSetType(sred->psubcomm, sred->subcommtype));
4729566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectMemory((PetscObject)pc, sizeof(PetscSubcomm)));
4738d9f7141SDave May         sred->subcomm = PetscSubcommChild(sred->psubcomm);
4748d9f7141SDave May       }
4758d9f7141SDave May     } else { /* query PC for DM, check communicators */
4768d9f7141SDave May       DM          dm, dm_coarse_partition          = NULL;
4778d9f7141SDave May       MPI_Comm    comm_fine, comm_coarse_partition = MPI_COMM_NULL;
4788d9f7141SDave May       PetscMPIInt csize_fine = 0, csize_coarse_partition = 0, cs[2], csg[2], cnt = 0;
4798d9f7141SDave May       PetscBool   isvalidsubcomm;
4808d9f7141SDave May 
4819566063dSJacob Faibussowitsch       PetscCall(PCGetDM(pc, &dm));
4828d9f7141SDave May       comm_fine = PetscObjectComm((PetscObject)dm);
4839566063dSJacob Faibussowitsch       PetscCall(DMGetCoarseDM(dm, &dm_coarse_partition));
4848d9f7141SDave May       if (dm_coarse_partition) { cnt = 1; }
4859566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &cnt, 1, MPI_INT, MPI_SUM, comm_fine));
48608401ef6SPierre Jolivet       PetscCheck(cnt != 0, comm_fine, PETSC_ERR_SUP, "Zero instances of a coarse DM were found");
4878d9f7141SDave May 
4889566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_size(comm_fine, &csize_fine));
4898d9f7141SDave May       if (dm_coarse_partition) {
4908d9f7141SDave May         comm_coarse_partition = PetscObjectComm((PetscObject)dm_coarse_partition);
4919566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(comm_coarse_partition, &csize_coarse_partition));
4928d9f7141SDave May       }
4938d9f7141SDave May 
4948d9f7141SDave May       cs[0] = csize_fine;
4958d9f7141SDave May       cs[1] = csize_coarse_partition;
4969566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allreduce(cs, csg, 2, MPI_INT, MPI_MAX, comm_fine));
49708401ef6SPierre Jolivet       PetscCheck(csg[0] != csg[1], comm_fine, PETSC_ERR_SUP, "Coarse DM uses the same size communicator as the parent DM attached to the PC");
4988d9f7141SDave May 
4999566063dSJacob Faibussowitsch       PetscCall(PCTelescopeTestValidSubcomm(comm_fine, comm_coarse_partition, &isvalidsubcomm));
50028b400f6SJacob Faibussowitsch       PetscCheck(isvalidsubcomm, comm_fine, PETSC_ERR_SUP, "Coarse DM communicator is not a sub-communicator of parentDM->comm");
5018d9f7141SDave May       sred->subcomm = comm_coarse_partition;
5028d9f7141SDave May     }
5038d9f7141SDave May   }
5048d9f7141SDave May   subcomm = sred->subcomm;
5058d9f7141SDave May 
5068d9f7141SDave May   /* internal KSP */
5078d9f7141SDave May   if (!pc->setupcalled) {
5088d9f7141SDave May     const char *prefix;
5098d9f7141SDave May 
51057f12427SDave May     if (PCTelescope_isActiveRank(sred)) {
5119566063dSJacob Faibussowitsch       PetscCall(KSPCreate(subcomm, &sred->ksp));
5129566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(sred->ksp, pc->erroriffailure));
5139566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)sred->ksp, (PetscObject)pc, 1));
5149566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)sred->ksp));
5159566063dSJacob Faibussowitsch       PetscCall(KSPSetType(sred->ksp, KSPPREONLY));
5169566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
5179566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(sred->ksp, prefix));
5189566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(sred->ksp, "telescope_"));
5198d9f7141SDave May     }
5201e07b27eSBarry Smith   }
5211e07b27eSBarry Smith 
5221e07b27eSBarry Smith   /* setup */
523*48a46eb9SPierre Jolivet   if (!pc->setupcalled && sred->pctelescope_setup_type) PetscCall(sred->pctelescope_setup_type(pc, sred));
5241e07b27eSBarry Smith   /* update */
5251e07b27eSBarry Smith   if (!pc->setupcalled) {
526*48a46eb9SPierre Jolivet     if (sred->pctelescope_matcreate_type) PetscCall(sred->pctelescope_matcreate_type(pc, sred, MAT_INITIAL_MATRIX, &sred->Bred));
5271baa6e33SBarry Smith     if (sred->pctelescope_matnullspacecreate_type) PetscCall(sred->pctelescope_matnullspacecreate_type(pc, sred, sred->Bred));
5281e07b27eSBarry Smith   } else {
529*48a46eb9SPierre Jolivet     if (sred->pctelescope_matcreate_type) PetscCall(sred->pctelescope_matcreate_type(pc, sred, MAT_REUSE_MATRIX, &sred->Bred));
5301e07b27eSBarry Smith   }
5311e07b27eSBarry Smith 
5321e07b27eSBarry Smith   /* common - no construction */
53357f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
5349566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(sred->ksp, sred->Bred, sred->Bred));
535*48a46eb9SPierre Jolivet     if (pc->setfromoptionscalled && !pc->setupcalled) PetscCall(KSPSetFromOptions(sred->ksp));
5361e07b27eSBarry Smith   }
5371e07b27eSBarry Smith   PetscFunctionReturn(0);
5381e07b27eSBarry Smith }
5391e07b27eSBarry Smith 
5409371c9d4SSatish Balay static PetscErrorCode PCApply_Telescope(PC pc, Vec x, Vec y) {
5411e07b27eSBarry Smith   PC_Telescope       sred = (PC_Telescope)pc->data;
5421e07b27eSBarry Smith   Vec                xtmp, xred, yred;
54313c30530SDave May   PetscInt           i, st, ed;
5441e07b27eSBarry Smith   VecScatter         scatter;
5451e07b27eSBarry Smith   PetscScalar       *array;
5461e07b27eSBarry Smith   const PetscScalar *x_array;
5471e07b27eSBarry Smith 
5481e07b27eSBarry Smith   PetscFunctionBegin;
5499566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
550bf00f589SPatrick Sanan 
5511e07b27eSBarry Smith   xtmp    = sred->xtmp;
5521e07b27eSBarry Smith   scatter = sred->scatter;
5531e07b27eSBarry Smith   xred    = sred->xred;
5541e07b27eSBarry Smith   yred    = sred->yred;
5551e07b27eSBarry Smith 
5561e07b27eSBarry Smith   /* pull in vector x->xtmp */
5579566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, x, xtmp, INSERT_VALUES, SCATTER_FORWARD));
5589566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, x, xtmp, INSERT_VALUES, SCATTER_FORWARD));
5591e07b27eSBarry Smith 
560bf00f589SPatrick Sanan   /* copy vector entries into xred */
5619566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xtmp, &x_array));
5621e07b27eSBarry Smith   if (xred) {
5631e07b27eSBarry Smith     PetscScalar *LA_xred;
5649566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(xred, &st, &ed));
5659566063dSJacob Faibussowitsch     PetscCall(VecGetArray(xred, &LA_xred));
5669371c9d4SSatish Balay     for (i = 0; i < ed - st; i++) { LA_xred[i] = x_array[i]; }
5679566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(xred, &LA_xred));
5681e07b27eSBarry Smith   }
5699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xtmp, &x_array));
5701e07b27eSBarry Smith   /* solve */
57157f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
5729566063dSJacob Faibussowitsch     PetscCall(KSPSolve(sred->ksp, xred, yred));
5739566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(sred->ksp, pc, yred));
5741e07b27eSBarry Smith   }
5751e07b27eSBarry Smith   /* return vector */
5769566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xtmp, &array));
5771e07b27eSBarry Smith   if (yred) {
5781e07b27eSBarry Smith     const PetscScalar *LA_yred;
5799566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(yred, &st, &ed));
5809566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(yred, &LA_yred));
5819371c9d4SSatish Balay     for (i = 0; i < ed - st; i++) { array[i] = LA_yred[i]; }
5829566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(yred, &LA_yred));
5831e07b27eSBarry Smith   }
5849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xtmp, &array));
5859566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, xtmp, y, INSERT_VALUES, SCATTER_REVERSE));
5869566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, xtmp, y, INSERT_VALUES, SCATTER_REVERSE));
5871e07b27eSBarry Smith   PetscFunctionReturn(0);
5881e07b27eSBarry Smith }
5891e07b27eSBarry Smith 
5909371c9d4SSatish Balay static PetscErrorCode PCApplyRichardson_Telescope(PC pc, Vec x, Vec y, Vec w, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt its, PetscBool zeroguess, PetscInt *outits, PCRichardsonConvergedReason *reason) {
591f650675bSDave May   PC_Telescope       sred = (PC_Telescope)pc->data;
592a1d91a28SDave May   Vec                xtmp, yred;
593f650675bSDave May   PetscInt           i, st, ed;
594f650675bSDave May   VecScatter         scatter;
595f650675bSDave May   const PetscScalar *x_array;
596f650675bSDave May   PetscBool          default_init_guess_value;
597f650675bSDave May 
598f650675bSDave May   PetscFunctionBegin;
599f650675bSDave May   xtmp    = sred->xtmp;
600f650675bSDave May   scatter = sred->scatter;
601f650675bSDave May   yred    = sred->yred;
602f650675bSDave May 
60308401ef6SPierre Jolivet   PetscCheck(its <= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "PCApplyRichardson_Telescope only supports max_it = 1");
604f650675bSDave May   *reason = (PCRichardsonConvergedReason)0;
605f650675bSDave May 
606f650675bSDave May   if (!zeroguess) {
6079566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "PCTelescope: Scattering y for non-zero initial guess\n"));
608f650675bSDave May     /* pull in vector y->xtmp */
6099566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scatter, y, xtmp, INSERT_VALUES, SCATTER_FORWARD));
6109566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scatter, y, xtmp, INSERT_VALUES, SCATTER_FORWARD));
611f650675bSDave May 
612bf00f589SPatrick Sanan     /* copy vector entries into xred */
6139566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(xtmp, &x_array));
614f650675bSDave May     if (yred) {
615f650675bSDave May       PetscScalar *LA_yred;
6169566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(yred, &st, &ed));
6179566063dSJacob Faibussowitsch       PetscCall(VecGetArray(yred, &LA_yred));
6189371c9d4SSatish Balay       for (i = 0; i < ed - st; i++) { LA_yred[i] = x_array[i]; }
6199566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(yred, &LA_yred));
620f650675bSDave May     }
6219566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(xtmp, &x_array));
622f650675bSDave May   }
623f650675bSDave May 
62457f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
6259566063dSJacob Faibussowitsch     PetscCall(KSPGetInitialGuessNonzero(sred->ksp, &default_init_guess_value));
6269566063dSJacob Faibussowitsch     if (!zeroguess) PetscCall(KSPSetInitialGuessNonzero(sred->ksp, PETSC_TRUE));
627f650675bSDave May   }
628f650675bSDave May 
6299566063dSJacob Faibussowitsch   PetscCall(PCApply_Telescope(pc, x, y));
630f650675bSDave May 
631*48a46eb9SPierre Jolivet   if (PCTelescope_isActiveRank(sred)) PetscCall(KSPSetInitialGuessNonzero(sred->ksp, default_init_guess_value));
632f650675bSDave May 
633f650675bSDave May   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
634f650675bSDave May   *outits = 1;
635f650675bSDave May   PetscFunctionReturn(0);
636f650675bSDave May }
637f650675bSDave May 
6389371c9d4SSatish Balay static PetscErrorCode PCReset_Telescope(PC pc) {
6391e07b27eSBarry Smith   PC_Telescope sred = (PC_Telescope)pc->data;
6401e07b27eSBarry Smith 
641362febeeSStefano Zampini   PetscFunctionBegin;
6429566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&sred->isin));
6439566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&sred->scatter));
6449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sred->xred));
6459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sred->yred));
6469566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sred->xtmp));
6479566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sred->Bred));
6489566063dSJacob Faibussowitsch   PetscCall(KSPReset(sred->ksp));
6491baa6e33SBarry Smith   if (sred->pctelescope_reset_type) PetscCall(sred->pctelescope_reset_type(pc));
6501e07b27eSBarry Smith   PetscFunctionReturn(0);
6511e07b27eSBarry Smith }
6521e07b27eSBarry Smith 
6539371c9d4SSatish Balay static PetscErrorCode PCDestroy_Telescope(PC pc) {
6541e07b27eSBarry Smith   PC_Telescope sred = (PC_Telescope)pc->data;
6551e07b27eSBarry Smith 
6561e07b27eSBarry Smith   PetscFunctionBegin;
6579566063dSJacob Faibussowitsch   PetscCall(PCReset_Telescope(pc));
6589566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&sred->ksp));
6599566063dSJacob Faibussowitsch   PetscCall(PetscSubcommDestroy(&sred->psubcomm));
6609566063dSJacob Faibussowitsch   PetscCall(PetscFree(sred->dm_ctx));
6612e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetKSP_C", NULL));
6622e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetSubcommType_C", NULL));
6632e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetSubcommType_C", NULL));
6642e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetReductionFactor_C", NULL));
6652e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetReductionFactor_C", NULL));
6662e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetIgnoreDM_C", NULL));
6672e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetIgnoreDM_C", NULL));
6682e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetIgnoreKSPComputeOperators_C", NULL));
6692e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetIgnoreKSPComputeOperators_C", NULL));
6702e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetDM_C", NULL));
6712e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetUseCoarseDM_C", NULL));
6722e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetUseCoarseDM_C", NULL));
6739566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
6741e07b27eSBarry Smith   PetscFunctionReturn(0);
6751e07b27eSBarry Smith }
6761e07b27eSBarry Smith 
6779371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_Telescope(PC pc, PetscOptionItems *PetscOptionsObject) {
6781e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
6791e07b27eSBarry Smith   MPI_Comm         comm;
6801e07b27eSBarry Smith   PetscMPIInt      size;
68148a10b22SPatrick Sanan   PetscBool        flg;
68248a10b22SPatrick Sanan   PetscSubcommType subcommtype;
6831e07b27eSBarry Smith 
6841e07b27eSBarry Smith   PetscFunctionBegin;
6859566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
6869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
687d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Telescope options");
6889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_telescope_subcomm_type", "Subcomm type (interlaced or contiguous)", "PCTelescopeSetSubcommType", PetscSubcommTypes, (PetscEnum)sred->subcommtype, (PetscEnum *)&subcommtype, &flg));
6891baa6e33SBarry Smith   if (flg) PetscCall(PCTelescopeSetSubcommType(pc, subcommtype));
6909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_telescope_reduction_factor", "Factor to reduce comm size by", "PCTelescopeSetReductionFactor", sred->redfactor, &sred->redfactor, NULL));
69108401ef6SPierre Jolivet   PetscCheck(sred->redfactor <= size, comm, PETSC_ERR_ARG_WRONG, "-pc_telescope_reduction_factor <= comm size");
6929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_telescope_ignore_dm", "Ignore any DM attached to the PC", "PCTelescopeSetIgnoreDM", sred->ignore_dm, &sred->ignore_dm, NULL));
6939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_telescope_ignore_kspcomputeoperators", "Ignore method used to compute A", "PCTelescopeSetIgnoreKSPComputeOperators", sred->ignore_kspcomputeoperators, &sred->ignore_kspcomputeoperators, NULL));
6949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_telescope_use_coarse_dm", "Define sub-communicator from the coarse DM", "PCTelescopeSetUseCoarseDM", sred->use_coarse_dm, &sred->use_coarse_dm, NULL));
695d0609cedSBarry Smith   PetscOptionsHeadEnd();
6961e07b27eSBarry Smith   PetscFunctionReturn(0);
6971e07b27eSBarry Smith }
6981e07b27eSBarry Smith 
6991e07b27eSBarry Smith /* PC simplementation specific API's */
7001e07b27eSBarry Smith 
7019371c9d4SSatish Balay static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc, KSP *ksp) {
7021e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
703bd49479cSSatish Balay   PetscFunctionBegin;
7041e07b27eSBarry Smith   if (ksp) *ksp = red->ksp;
705bd49479cSSatish Balay   PetscFunctionReturn(0);
7061e07b27eSBarry Smith }
7071e07b27eSBarry Smith 
7089371c9d4SSatish Balay static PetscErrorCode PCTelescopeGetSubcommType_Telescope(PC pc, PetscSubcommType *subcommtype) {
70948a10b22SPatrick Sanan   PC_Telescope red = (PC_Telescope)pc->data;
71048a10b22SPatrick Sanan   PetscFunctionBegin;
71148a10b22SPatrick Sanan   if (subcommtype) *subcommtype = red->subcommtype;
71248a10b22SPatrick Sanan   PetscFunctionReturn(0);
71348a10b22SPatrick Sanan }
71448a10b22SPatrick Sanan 
7159371c9d4SSatish Balay static PetscErrorCode PCTelescopeSetSubcommType_Telescope(PC pc, PetscSubcommType subcommtype) {
71648a10b22SPatrick Sanan   PC_Telescope red = (PC_Telescope)pc->data;
71748a10b22SPatrick Sanan 
71848a10b22SPatrick Sanan   PetscFunctionBegin;
71928b400f6SJacob Faibussowitsch   PetscCheck(!pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "You cannot change the subcommunicator type for PCTelescope after it has been set up.");
72048a10b22SPatrick Sanan   red->subcommtype = subcommtype;
72148a10b22SPatrick Sanan   PetscFunctionReturn(0);
72248a10b22SPatrick Sanan }
72348a10b22SPatrick Sanan 
7249371c9d4SSatish Balay static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc, PetscInt *fact) {
7251e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
726bd49479cSSatish Balay   PetscFunctionBegin;
7271e07b27eSBarry Smith   if (fact) *fact = red->redfactor;
728bd49479cSSatish Balay   PetscFunctionReturn(0);
7291e07b27eSBarry Smith }
7301e07b27eSBarry Smith 
7319371c9d4SSatish Balay static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc, PetscInt fact) {
7321e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
7331e07b27eSBarry Smith   PetscMPIInt  size;
7341e07b27eSBarry Smith 
735bd49479cSSatish Balay   PetscFunctionBegin;
7369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
73763a3b9bcSJacob Faibussowitsch   PetscCheck(fact > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Reduction factor of telescoping PC %" PetscInt_FMT " must be positive", fact);
73863a3b9bcSJacob Faibussowitsch   PetscCheck(fact <= size, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Reduction factor of telescoping PC %" PetscInt_FMT " must be <= comm.size", fact);
7391e07b27eSBarry Smith   red->redfactor = fact;
740bd49479cSSatish Balay   PetscFunctionReturn(0);
7411e07b27eSBarry Smith }
7421e07b27eSBarry Smith 
7439371c9d4SSatish Balay static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc, PetscBool *v) {
7441e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
745bd49479cSSatish Balay   PetscFunctionBegin;
7461e07b27eSBarry Smith   if (v) *v = red->ignore_dm;
747bd49479cSSatish Balay   PetscFunctionReturn(0);
7481e07b27eSBarry Smith }
74948a10b22SPatrick Sanan 
7509371c9d4SSatish Balay static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc, PetscBool v) {
7511e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
752bd49479cSSatish Balay   PetscFunctionBegin;
7531e07b27eSBarry Smith   red->ignore_dm = v;
754bd49479cSSatish Balay   PetscFunctionReturn(0);
7551e07b27eSBarry Smith }
7561e07b27eSBarry Smith 
7579371c9d4SSatish Balay static PetscErrorCode PCTelescopeGetUseCoarseDM_Telescope(PC pc, PetscBool *v) {
7588d9f7141SDave May   PC_Telescope red = (PC_Telescope)pc->data;
7598d9f7141SDave May   PetscFunctionBegin;
7608d9f7141SDave May   if (v) *v = red->use_coarse_dm;
7618d9f7141SDave May   PetscFunctionReturn(0);
7628d9f7141SDave May }
7638d9f7141SDave May 
7649371c9d4SSatish Balay static PetscErrorCode PCTelescopeSetUseCoarseDM_Telescope(PC pc, PetscBool v) {
7658d9f7141SDave May   PC_Telescope red = (PC_Telescope)pc->data;
7668d9f7141SDave May   PetscFunctionBegin;
7678d9f7141SDave May   red->use_coarse_dm = v;
7688d9f7141SDave May   PetscFunctionReturn(0);
7698d9f7141SDave May }
7708d9f7141SDave May 
7719371c9d4SSatish Balay static PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc, PetscBool *v) {
7720ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
7730ae7c45bSDave May   PetscFunctionBegin;
7740ae7c45bSDave May   if (v) *v = red->ignore_kspcomputeoperators;
7750ae7c45bSDave May   PetscFunctionReturn(0);
7760ae7c45bSDave May }
77748a10b22SPatrick Sanan 
7789371c9d4SSatish Balay static PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc, PetscBool v) {
7790ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
7800ae7c45bSDave May   PetscFunctionBegin;
7810ae7c45bSDave May   red->ignore_kspcomputeoperators = v;
7820ae7c45bSDave May   PetscFunctionReturn(0);
7830ae7c45bSDave May }
7840ae7c45bSDave May 
7859371c9d4SSatish Balay static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc, DM *dm) {
7861e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
787bd49479cSSatish Balay   PetscFunctionBegin;
7881e07b27eSBarry Smith   *dm = private_PCTelescopeGetSubDM(red);
789bd49479cSSatish Balay   PetscFunctionReturn(0);
7901e07b27eSBarry Smith }
7911e07b27eSBarry Smith 
7921e07b27eSBarry Smith /*@
7931e07b27eSBarry Smith  PCTelescopeGetKSP - Gets the KSP created by the telescoping PC.
7941e07b27eSBarry Smith 
7951e07b27eSBarry Smith  Not Collective
7961e07b27eSBarry Smith 
7971e07b27eSBarry Smith  Input Parameter:
7981e07b27eSBarry Smith .  pc - the preconditioner context
7991e07b27eSBarry Smith 
8001e07b27eSBarry Smith  Output Parameter:
8011e07b27eSBarry Smith .  subksp - the KSP defined the smaller set of processes
8021e07b27eSBarry Smith 
8031e07b27eSBarry Smith  Level: advanced
8041e07b27eSBarry Smith 
8051e07b27eSBarry Smith @*/
8069371c9d4SSatish Balay PetscErrorCode PCTelescopeGetKSP(PC pc, KSP *subksp) {
807bd49479cSSatish Balay   PetscFunctionBegin;
808cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetKSP_C", (PC, KSP *), (pc, subksp));
809bd49479cSSatish Balay   PetscFunctionReturn(0);
8101e07b27eSBarry Smith }
8111e07b27eSBarry Smith 
8121e07b27eSBarry Smith /*@
8131e07b27eSBarry Smith  PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by.
8141e07b27eSBarry Smith 
8151e07b27eSBarry Smith  Not Collective
8161e07b27eSBarry Smith 
8171e07b27eSBarry Smith  Input Parameter:
8181e07b27eSBarry Smith .  pc - the preconditioner context
8191e07b27eSBarry Smith 
8201e07b27eSBarry Smith  Output Parameter:
8211e07b27eSBarry Smith .  fact - the reduction factor
8221e07b27eSBarry Smith 
8231e07b27eSBarry Smith  Level: advanced
8241e07b27eSBarry Smith 
8251e07b27eSBarry Smith @*/
8269371c9d4SSatish Balay PetscErrorCode PCTelescopeGetReductionFactor(PC pc, PetscInt *fact) {
827bd49479cSSatish Balay   PetscFunctionBegin;
828cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetReductionFactor_C", (PC, PetscInt *), (pc, fact));
829bd49479cSSatish Balay   PetscFunctionReturn(0);
8301e07b27eSBarry Smith }
8311e07b27eSBarry Smith 
8321e07b27eSBarry Smith /*@
8331e07b27eSBarry Smith  PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by.
8341e07b27eSBarry Smith 
8351e07b27eSBarry Smith  Not Collective
8361e07b27eSBarry Smith 
8371e07b27eSBarry Smith  Input Parameter:
8381e07b27eSBarry Smith .  pc - the preconditioner context
8391e07b27eSBarry Smith 
8401e07b27eSBarry Smith  Output Parameter:
8411e07b27eSBarry Smith .  fact - the reduction factor
8421e07b27eSBarry Smith 
8431e07b27eSBarry Smith  Level: advanced
8441e07b27eSBarry Smith 
8451e07b27eSBarry Smith @*/
8469371c9d4SSatish Balay PetscErrorCode PCTelescopeSetReductionFactor(PC pc, PetscInt fact) {
847bd49479cSSatish Balay   PetscFunctionBegin;
848cac4c232SBarry Smith   PetscTryMethod(pc, "PCTelescopeSetReductionFactor_C", (PC, PetscInt), (pc, fact));
849bd49479cSSatish Balay   PetscFunctionReturn(0);
8501e07b27eSBarry Smith }
8511e07b27eSBarry Smith 
8521e07b27eSBarry Smith /*@
8531e07b27eSBarry Smith  PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used.
8541e07b27eSBarry Smith 
8551e07b27eSBarry Smith  Not Collective
8561e07b27eSBarry Smith 
8571e07b27eSBarry Smith  Input Parameter:
8581e07b27eSBarry Smith .  pc - the preconditioner context
8591e07b27eSBarry Smith 
8601e07b27eSBarry Smith  Output Parameter:
8611e07b27eSBarry Smith .  v - the flag
8621e07b27eSBarry Smith 
8631e07b27eSBarry Smith  Level: advanced
8641e07b27eSBarry Smith 
8651e07b27eSBarry Smith @*/
8669371c9d4SSatish Balay PetscErrorCode PCTelescopeGetIgnoreDM(PC pc, PetscBool *v) {
867bd49479cSSatish Balay   PetscFunctionBegin;
868cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetIgnoreDM_C", (PC, PetscBool *), (pc, v));
869bd49479cSSatish Balay   PetscFunctionReturn(0);
8701e07b27eSBarry Smith }
8711e07b27eSBarry Smith 
8721e07b27eSBarry Smith /*@
8731e07b27eSBarry Smith  PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC.
8741e07b27eSBarry Smith 
8751e07b27eSBarry Smith  Not Collective
8761e07b27eSBarry Smith 
8771e07b27eSBarry Smith  Input Parameter:
8781e07b27eSBarry Smith .  pc - the preconditioner context
8791e07b27eSBarry Smith 
8801e07b27eSBarry Smith  Output Parameter:
8811e07b27eSBarry Smith .  v - Use PETSC_TRUE to ignore any DM
8821e07b27eSBarry Smith 
8831e07b27eSBarry Smith  Level: advanced
8841e07b27eSBarry Smith 
8851e07b27eSBarry Smith @*/
8869371c9d4SSatish Balay PetscErrorCode PCTelescopeSetIgnoreDM(PC pc, PetscBool v) {
887bd49479cSSatish Balay   PetscFunctionBegin;
888cac4c232SBarry Smith   PetscTryMethod(pc, "PCTelescopeSetIgnoreDM_C", (PC, PetscBool), (pc, v));
889bd49479cSSatish Balay   PetscFunctionReturn(0);
8901e07b27eSBarry Smith }
8911e07b27eSBarry Smith 
8921e07b27eSBarry Smith /*@
8938d9f7141SDave May  PCTelescopeGetUseCoarseDM - Get the flag indicating if the coarse DM attached to DM associated with the PC will be used.
8948d9f7141SDave May 
8958d9f7141SDave May  Not Collective
8968d9f7141SDave May 
8978d9f7141SDave May  Input Parameter:
8988d9f7141SDave May .  pc - the preconditioner context
8998d9f7141SDave May 
9008d9f7141SDave May  Output Parameter:
9018d9f7141SDave May .  v - the flag
9028d9f7141SDave May 
9038d9f7141SDave May  Level: advanced
9048d9f7141SDave May 
9058d9f7141SDave May @*/
9069371c9d4SSatish Balay PetscErrorCode PCTelescopeGetUseCoarseDM(PC pc, PetscBool *v) {
9078d9f7141SDave May   PetscFunctionBegin;
908cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetUseCoarseDM_C", (PC, PetscBool *), (pc, v));
9098d9f7141SDave May   PetscFunctionReturn(0);
9108d9f7141SDave May }
9118d9f7141SDave May 
9128d9f7141SDave May /*@
9138d9f7141SDave May  PCTelescopeSetUseCoarseDM - Set a flag to query the DM attached to the PC if it also has a coarse DM
9148d9f7141SDave May 
9158d9f7141SDave May  Not Collective
9168d9f7141SDave May 
9178d9f7141SDave May  Input Parameter:
9188d9f7141SDave May .  pc - the preconditioner context
9198d9f7141SDave May 
9208d9f7141SDave May  Output Parameter:
921aaa7a805SDave May .  v - Use PETSC_FALSE to ignore any coarse DM
9228d9f7141SDave May 
9238d9f7141SDave May  Notes:
9248d9f7141SDave May  When you have specified to use a coarse DM, the communicator used to create the sub-KSP within PCTelescope
9258d9f7141SDave May  will be that of the coarse DM. Hence the flags -pc_telescope_reduction_factor and
9268d9f7141SDave May  -pc_telescope_subcomm_type will no longer have any meaning.
9278d9f7141SDave May  It is required that the communicator associated with the parent (fine) and the coarse DM are of different sizes.
9288d9f7141SDave May  An error will occur of the size of the communicator associated with the coarse DM
9298d9f7141SDave May  is the same as that of the parent DM.
9308d9f7141SDave May  Furthermore, it is required that the communicator on the coarse DM is a sub-communicator of the parent.
9318d9f7141SDave May  This will be checked at the time the preconditioner is setup and an error will occur if
9328d9f7141SDave May  the coarse DM does not define a sub-communicator of that used by the parent DM.
9338d9f7141SDave May 
9348d9f7141SDave May  The particular Telescope setup invoked when using a coarse DM is agnostic with respect to the type of
9358d9f7141SDave May  the DM used (e.g. it supports DMSHELL, DMPLEX, etc).
9368d9f7141SDave May 
9378d9f7141SDave May  Support is currently only provided for the case when you are using KSPSetComputeOperators()
9388d9f7141SDave May 
9398d9f7141SDave May  The user is required to compose a function with the parent DM to facilitate the transfer of fields (Vec) between the different decompositions defined by the fine and coarse DMs.
9408d9f7141SDave May  In the user code, this is achieved via
9418d9f7141SDave May .vb
9428d9f7141SDave May    {
9438d9f7141SDave May      DM dm_fine;
9448d9f7141SDave May      PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeFieldScatter",your_field_scatter_method);
9458d9f7141SDave May    }
9468d9f7141SDave May .ve
9478d9f7141SDave May  The signature of the user provided field scatter method is
9488d9f7141SDave May .vb
9498d9f7141SDave May    PetscErrorCode your_field_scatter_method(DM dm_fine,Vec x_fine,ScatterMode mode,DM dm_coarse,Vec x_coarse);
9508d9f7141SDave May .ve
9518d9f7141SDave May  The user must provide support for both mode = SCATTER_FORWARD and mode = SCATTER_REVERSE.
9528d9f7141SDave May  SCATTER_FORWARD implies the direction of transfer is from the parent (fine) DM to the coarse DM.
9538d9f7141SDave May 
9548d9f7141SDave May  Optionally, the user may also compose a function with the parent DM to facilitate the transfer
9558d9f7141SDave May  of state variables between the fine and coarse DMs.
9568d9f7141SDave May  In the context of a finite element discretization, an example state variable might be
9578d9f7141SDave May  values associated with quadrature points within each element.
9588d9f7141SDave May  A user provided state scatter method is composed via
9598d9f7141SDave May .vb
9608d9f7141SDave May    {
9618d9f7141SDave May      DM dm_fine;
9628d9f7141SDave May      PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeStateScatter",your_state_scatter_method);
9638d9f7141SDave May    }
9648d9f7141SDave May .ve
9658d9f7141SDave May  The signature of the user provided state scatter method is
9668d9f7141SDave May .vb
9678d9f7141SDave May    PetscErrorCode your_state_scatter_method(DM dm_fine,ScatterMode mode,DM dm_coarse);
9688d9f7141SDave May .ve
9698d9f7141SDave May  SCATTER_FORWARD implies the direction of transfer is from the fine DM to the coarse DM.
9708d9f7141SDave May  The user is only required to support mode = SCATTER_FORWARD.
9718d9f7141SDave May  No assumption is made about the data type of the state variables.
9728d9f7141SDave May  These must be managed by the user and must be accessible from the DM.
9738d9f7141SDave May 
9748d9f7141SDave May  Care must be taken in defining the user context passed to KSPSetComputeOperators() which is to be
9758d9f7141SDave May  associated with the sub-KSP residing within PCTelescope.
9768d9f7141SDave May  In general, PCTelescope assumes that the context on the fine and coarse DM used with
9778d9f7141SDave May  KSPSetComputeOperators() should be "similar" in type or origin.
9788d9f7141SDave May  Specifically the following rules are used to infer what context on the sub-KSP should be.
9798d9f7141SDave May 
9808d9f7141SDave May  First the contexts from the KSP and the fine and coarse DMs are retrieved.
9818d9f7141SDave May  Note that the special case of a DMSHELL context is queried.
9828d9f7141SDave May 
9838d9f7141SDave May .vb
9848d9f7141SDave May    DMKSPGetComputeOperators(dm_fine,&dmfine_kspfunc,&dmfine_kspctx);
9858d9f7141SDave May    DMGetApplicationContext(dm_fine,&dmfine_appctx);
9868d9f7141SDave May    DMShellGetContext(dm_fine,&dmfine_shellctx);
9878d9f7141SDave May 
9888d9f7141SDave May    DMGetApplicationContext(dm_coarse,&dmcoarse_appctx);
9898d9f7141SDave May    DMShellGetContext(dm_coarse,&dmcoarse_shellctx);
9908d9f7141SDave May .ve
9918d9f7141SDave May 
9928d9f7141SDave May  The following rules are then enforced:
9938d9f7141SDave May 
9948d9f7141SDave May  1. If dmfine_kspctx = NULL, then we provide a NULL pointer as the context for the sub-KSP:
9958d9f7141SDave May  KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,NULL);
9968d9f7141SDave May 
9978d9f7141SDave May  2. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_appctx,
9988d9f7141SDave May  check that dmcoarse_appctx is also non-NULL. If this is true, then:
9998d9f7141SDave May  KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,dmcoarse_appctx);
10008d9f7141SDave May 
10018d9f7141SDave May  3. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_shellctx,
10028d9f7141SDave May  check that dmcoarse_shellctx is also non-NULL. If this is true, then:
10038d9f7141SDave May  KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,dmcoarse_shellctx);
10048d9f7141SDave May 
10058d9f7141SDave May  If neither of the above three tests passed, then PCTelescope cannot safely determine what
10068d9f7141SDave May  context should be provided to KSPSetComputeOperators() for use with the sub-KSP.
10078d9f7141SDave May  In this case, an additional mechanism is provided via a composed function which will return
10088d9f7141SDave May  the actual context to be used. To use this feature you must compose the "getter" function
10098d9f7141SDave May  with the coarse DM, e.g.
10108d9f7141SDave May .vb
10118d9f7141SDave May    {
10128d9f7141SDave May      DM dm_coarse;
10138d9f7141SDave May      PetscObjectCompose((PetscObject)dm_coarse,"PCTelescopeGetCoarseDMKSPContext",your_coarse_context_getter);
10148d9f7141SDave May    }
10158d9f7141SDave May .ve
10168d9f7141SDave May  The signature of the user provided method is
10178d9f7141SDave May .vb
10188d9f7141SDave May    PetscErrorCode your_coarse_context_getter(DM dm_coarse,void **your_kspcontext);
10198d9f7141SDave May .ve
10208d9f7141SDave May 
10218d9f7141SDave May  Level: advanced
10228d9f7141SDave May 
10238d9f7141SDave May @*/
10249371c9d4SSatish Balay PetscErrorCode PCTelescopeSetUseCoarseDM(PC pc, PetscBool v) {
10258d9f7141SDave May   PetscFunctionBegin;
1026cac4c232SBarry Smith   PetscTryMethod(pc, "PCTelescopeSetUseCoarseDM_C", (PC, PetscBool), (pc, v));
10278d9f7141SDave May   PetscFunctionReturn(0);
10288d9f7141SDave May }
10298d9f7141SDave May 
10308d9f7141SDave May /*@
10310ae7c45bSDave May  PCTelescopeGetIgnoreKSPComputeOperators - Get the flag indicating if KSPComputeOperators will be used.
10320ae7c45bSDave May 
10330ae7c45bSDave May  Not Collective
10340ae7c45bSDave May 
10350ae7c45bSDave May  Input Parameter:
10360ae7c45bSDave May .  pc - the preconditioner context
10370ae7c45bSDave May 
10380ae7c45bSDave May  Output Parameter:
10390ae7c45bSDave May .  v - the flag
10400ae7c45bSDave May 
10410ae7c45bSDave May  Level: advanced
10420ae7c45bSDave May 
10430ae7c45bSDave May @*/
10449371c9d4SSatish Balay PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC pc, PetscBool *v) {
10450ae7c45bSDave May   PetscFunctionBegin;
1046cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetIgnoreKSPComputeOperators_C", (PC, PetscBool *), (pc, v));
10470ae7c45bSDave May   PetscFunctionReturn(0);
10480ae7c45bSDave May }
10490ae7c45bSDave May 
10500ae7c45bSDave May /*@
10510ae7c45bSDave May  PCTelescopeSetIgnoreKSPComputeOperators - Set a flag to ignore KSPComputeOperators.
10520ae7c45bSDave May 
10530ae7c45bSDave May  Not Collective
10540ae7c45bSDave May 
10550ae7c45bSDave May  Input Parameter:
10560ae7c45bSDave May .  pc - the preconditioner context
10570ae7c45bSDave May 
10580ae7c45bSDave May  Output Parameter:
1059a954d8f4SDave May .  v - Use PETSC_TRUE to ignore the method (if defined) set via KSPSetComputeOperators on pc
10600ae7c45bSDave May 
10610ae7c45bSDave May  Level: advanced
10620ae7c45bSDave May 
10630ae7c45bSDave May @*/
10649371c9d4SSatish Balay PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC pc, PetscBool v) {
10650ae7c45bSDave May   PetscFunctionBegin;
1066cac4c232SBarry Smith   PetscTryMethod(pc, "PCTelescopeSetIgnoreKSPComputeOperators_C", (PC, PetscBool), (pc, v));
10670ae7c45bSDave May   PetscFunctionReturn(0);
10680ae7c45bSDave May }
10690ae7c45bSDave May 
10700ae7c45bSDave May /*@
10711e07b27eSBarry Smith  PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP.
10721e07b27eSBarry Smith 
10731e07b27eSBarry Smith  Not Collective
10741e07b27eSBarry Smith 
10751e07b27eSBarry Smith  Input Parameter:
10761e07b27eSBarry Smith .  pc - the preconditioner context
10771e07b27eSBarry Smith 
10781e07b27eSBarry Smith  Output Parameter:
10791e07b27eSBarry Smith .  subdm - The re-partitioned DM
10801e07b27eSBarry Smith 
10811e07b27eSBarry Smith  Level: advanced
10821e07b27eSBarry Smith 
10831e07b27eSBarry Smith @*/
10849371c9d4SSatish Balay PetscErrorCode PCTelescopeGetDM(PC pc, DM *subdm) {
1085bd49479cSSatish Balay   PetscFunctionBegin;
1086cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetDM_C", (PC, DM *), (pc, subdm));
1087bd49479cSSatish Balay   PetscFunctionReturn(0);
10881e07b27eSBarry Smith }
10891e07b27eSBarry Smith 
109048a10b22SPatrick Sanan /*@
109148a10b22SPatrick Sanan  PCTelescopeSetSubcommType - set subcommunicator type (interlaced or contiguous)
109248a10b22SPatrick Sanan 
109348a10b22SPatrick Sanan  Logically Collective
109448a10b22SPatrick Sanan 
1095d8d19677SJose E. Roman  Input Parameters:
10961dae98e4SBarry Smith +  pc - the preconditioner context
10971dae98e4SBarry Smith -  subcommtype - the subcommunicator type (see PetscSubcommType)
109848a10b22SPatrick Sanan 
109948a10b22SPatrick Sanan  Level: advanced
110048a10b22SPatrick Sanan 
1101db781477SPatrick Sanan .seealso: `PetscSubcommType`, `PetscSubcomm`, `PCTELESCOPE`
110248a10b22SPatrick Sanan @*/
11039371c9d4SSatish Balay PetscErrorCode PCTelescopeSetSubcommType(PC pc, PetscSubcommType subcommtype) {
110448a10b22SPatrick Sanan   PetscFunctionBegin;
1105cac4c232SBarry Smith   PetscTryMethod(pc, "PCTelescopeSetSubcommType_C", (PC, PetscSubcommType), (pc, subcommtype));
110648a10b22SPatrick Sanan   PetscFunctionReturn(0);
110748a10b22SPatrick Sanan }
110848a10b22SPatrick Sanan 
110948a10b22SPatrick Sanan /*@
111048a10b22SPatrick Sanan  PCTelescopeGetSubcommType - Get the subcommunicator type (interlaced or contiguous)
111148a10b22SPatrick Sanan 
111248a10b22SPatrick Sanan  Not Collective
111348a10b22SPatrick Sanan 
111448a10b22SPatrick Sanan  Input Parameter:
111548a10b22SPatrick Sanan .  pc - the preconditioner context
111648a10b22SPatrick Sanan 
111748a10b22SPatrick Sanan  Output Parameter:
111848a10b22SPatrick Sanan .  subcommtype - the subcommunicator type (see PetscSubcommType)
111948a10b22SPatrick Sanan 
112048a10b22SPatrick Sanan  Level: advanced
112148a10b22SPatrick Sanan 
1122db781477SPatrick Sanan .seealso: `PetscSubcomm`, `PetscSubcommType`, `PCTELESCOPE`
112348a10b22SPatrick Sanan @*/
11249371c9d4SSatish Balay PetscErrorCode PCTelescopeGetSubcommType(PC pc, PetscSubcommType *subcommtype) {
112548a10b22SPatrick Sanan   PetscFunctionBegin;
1126cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetSubcommType_C", (PC, PetscSubcommType *), (pc, subcommtype));
112748a10b22SPatrick Sanan   PetscFunctionReturn(0);
112848a10b22SPatrick Sanan }
112948a10b22SPatrick Sanan 
11301e07b27eSBarry Smith /* -------------------------------------------------------------------------------------*/
11311e07b27eSBarry Smith /*MC
113200fea0ebSDave May    PCTELESCOPE - Runs a KSP solver on a sub-communicator. MPI ranks not in the sub-communicator are idle during the solve.
11331e07b27eSBarry Smith 
11341e07b27eSBarry Smith    Options Database:
113500fea0ebSDave May +  -pc_telescope_reduction_factor <r> - factor to reduce the communicator size by. e.g. with 64 MPI ranks and r=4, the new sub-communicator will have 64/4 = 16 ranks.
113600fea0ebSDave May .  -pc_telescope_ignore_dm  - flag to indicate whether an attached DM should be ignored.
113700fea0ebSDave May .  -pc_telescope_subcomm_type <interlaced,contiguous> - defines the selection of MPI ranks on the sub-communicator. see PetscSubcomm for more information.
113800fea0ebSDave May .  -pc_telescope_ignore_kspcomputeoperators - flag to indicate whether KSPSetComputeOperators should be used on the sub-KSP.
113900fea0ebSDave May -  -pc_telescope_use_coarse_dm - flag to indicate whether the coarse DM should be used to define the sub-communicator.
11401e07b27eSBarry Smith 
11411e07b27eSBarry Smith    Level: advanced
11421e07b27eSBarry Smith 
11431e07b27eSBarry Smith    Notes:
114400fea0ebSDave May    Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation
114500fea0ebSDave May    creates a child sub-communicator (c') containing fewer MPI ranks than the original parent preconditioner (PC).
11466fc41876SBarry Smith    The preconditioner is deemed telescopic as it only calls KSPSolve() on a single
11478439623fSPatrick Sanan    sub-communicator, in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators.
114800fea0ebSDave May    This means there will be MPI ranks which will be idle during the application of this preconditioner.
114900fea0ebSDave May    Additionally, in comparison with PCREDUNDANT, PCTELESCOPE can utilize an attached DM.
11506fc41876SBarry Smith 
115100fea0ebSDave May    The default type of the sub KSP (the KSP defined on c') is PREONLY.
115200fea0ebSDave May 
115300fea0ebSDave May    There are three setup mechanisms for PCTelescope. Features support by each type are described below.
115400fea0ebSDave May    In the following, we will refer to the operators B and B', these are the Bmat provided to the KSP on the
115500fea0ebSDave May    communicators c and c' respectively.
115600fea0ebSDave May 
115700fea0ebSDave May    [1] Default setup
115800fea0ebSDave May    The sub-communicator c' is created via PetscSubcommCreate().
1159a5b23f4aSJose E. Roman    Explicitly defined nullspace and near nullspace vectors will be propagated from B to B'.
116000fea0ebSDave May    Currently there is no support define nullspaces via a user supplied method (e.g. as passed to MatNullSpaceSetFunction()).
116100fea0ebSDave May    No support is provided for KSPSetComputeOperators().
116200fea0ebSDave May    Currently there is no support for the flag -pc_use_amat.
116300fea0ebSDave May 
116400fea0ebSDave May    [2] DM aware setup
116500fea0ebSDave May    If a DM is attached to the PC, it is re-partitioned on the sub-communicator c'.
116600fea0ebSDave May    c' is created via PetscSubcommCreate().
11671e07b27eSBarry Smith    Both the Bmat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned DM.
11681e07b27eSBarry Smith    Currently only support for re-partitioning a DMDA is provided.
116900fea0ebSDave May    Any explicitly defined nullspace or near nullspace vectors attached to the original Bmat operator (B) are extracted, re-partitioned and set on the re-partitioned Bmat operator (B').
117000fea0ebSDave May    Currently there is no support define nullspaces via a user supplied method (e.g. as passed to MatNullSpaceSetFunction()).
117100fea0ebSDave May    Support is provided for KSPSetComputeOperators(). The user provided function and context is propagated to the sub KSP.
117200fea0ebSDave May    This is fragile since the user must ensure that their user context is valid for use on c'.
117300fea0ebSDave May    Currently there is no support for the flag -pc_use_amat.
11741e07b27eSBarry Smith 
117500fea0ebSDave May    [3] Coarse DM setup
117600fea0ebSDave May    If a DM (dmfine) is attached to the PC, dmfine is queried for a "coarse" DM (call this dmcoarse) via DMGetCoarseDM().
117700fea0ebSDave May    PCTELESCOPE will interpret the coarse DM as being defined on a sub-communicator of c.
117800fea0ebSDave May    The communicator associated with dmcoarse will define the c' to be used within PCTELESCOPE.
117900fea0ebSDave May    PCTELESCOPE will check that c' is in fact a sub-communicator of c. If it is not, an error will be reported.
118000fea0ebSDave May    The intention of this setup type is that PCTELESCOPE will use an existing (e.g. user defined) communicator hierarchy, say as would be
118100fea0ebSDave May    available with using multi-grid on unstructured meshes.
118200fea0ebSDave May    This setup will not use the command line options -pc_telescope_reduction_factor or -pc_telescope_subcomm_type.
118300fea0ebSDave May    Any explicitly defined nullspace or near nullspace vectors attached to the original Bmat operator (B) are extracted, scattered into the correct ordering consistent with dmcoarse and set on B'.
118400fea0ebSDave May    Currently there is no support define nullspaces via a user supplied method (e.g. as passed to MatNullSpaceSetFunction()).
118500fea0ebSDave May    There is no general method to permute field orderings, hence only KSPSetComputeOperators() is supported.
118600fea0ebSDave May    The user must use PetscObjectComposeFunction() with dmfine to define the method to scatter fields from dmfine to dmcoarse.
1187a5b23f4aSJose E. Roman    Propagation of the user context for KSPSetComputeOperators() on the sub KSP is attempted by querying the DM contexts associated with dmfine and dmcoarse. Alternatively, the user may use PetscObjectComposeFunction() with dmcoarse to define a method which will return the appropriate user context for KSPSetComputeOperators().
118800fea0ebSDave May    Currently there is no support for the flag -pc_use_amat.
118900fea0ebSDave May    This setup can be invoked by the option -pc_telescope_use_coarse_dm or by calling PCTelescopeSetUseCoarseDM(pc,PETSC_TRUE);
119000fea0ebSDave May    Further information about the user-provided methods required by this setup type are described here PCTelescopeSetUseCoarseDM().
11916fc41876SBarry Smith 
11926fc41876SBarry Smith    Developer Notes:
11936fc41876SBarry Smith    During PCSetup, the B operator is scattered onto c'.
11946fc41876SBarry Smith    Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c').
11958439623fSPatrick Sanan    Then, KSPSolve() is executed on the c' communicator.
11966fc41876SBarry Smith 
11976fc41876SBarry Smith    The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED
1198a04a6428SPatrick Sanan    creation routine by default (this can be changed with -pc_telescope_subcomm_type). We run the sub KSP on only the ranks within the communicator which have a color equal to zero.
11996fc41876SBarry Smith 
1200005d9f20SPatrick Sanan    The telescoping preconditioner is aware of nullspaces and near nullspaces which are attached to the B operator.
12018439623fSPatrick Sanan    In the case where B has a (near) nullspace attached, the (near) nullspace vectors are extracted from B and mapped into
1202005d9f20SPatrick Sanan    a new (near) nullspace, defined on the sub-communicator, which is attached to B' (the B operator which was scattered to c')
12036fc41876SBarry Smith 
120400fea0ebSDave May    The telescoping preconditioner can re-partition an attached DM if it is a DMDA (2D or 3D -
12054c500f23SPierre Jolivet    support for 1D DMDAs is not provided). If a DMDA is found, a topologically equivalent DMDA is created on c'
120600fea0ebSDave May    and this new DM is attached the sub KSP. The design of telescope is such that it should be possible to extend support
12078439623fSPatrick Sanan    for re-partitioning other to DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs.
120800fea0ebSDave May    Alternatively, user-provided re-partitioned DMs can be used via -pc_telescope_use_coarse_dm.
12096fc41876SBarry Smith 
121000fea0ebSDave May    With the default setup mode, B' is defined by fusing rows (in order) associated with MPI ranks common to c and c'.
12116fc41876SBarry Smith 
12128439623fSPatrick Sanan    When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permutation of B
12137dae84e0SHong Zhang    into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatCreateSubMatrices(), (iii) fusing the
12146fc41876SBarry Smith    locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat()
12156fc41876SBarry Smith 
12168439623fSPatrick Sanan    Limitations/improvements include the following.
12178439623fSPatrick Sanan    VecPlaceArray() could be used within PCApply() to improve efficiency and reduce memory usage.
121800fea0ebSDave May    A unified mechanism to query for user contexts as required by KSPSetComputeOperators() and MatNullSpaceSetFunction().
12196fc41876SBarry Smith 
12200705ae88SPierre Jolivet    The symmetric permutation used when a DMDA is encountered is performed via explicitly assembling a permutation matrix P,
12218439623fSPatrick Sanan    and performing P^T.A.P. Possibly it might be more efficient to use MatPermute(). We opted to use P^T.A.P as it appears
12220705ae88SPierre Jolivet    VecPermute() does not support the use case required here. By computing P, one can permute both the operator and RHS in a
12236fc41876SBarry Smith    consistent manner.
12246fc41876SBarry Smith 
122500fea0ebSDave May    Mapping of vectors (default setup mode) is performed in the following way.
122600fea0ebSDave May    Suppose the parent communicator size was 4, and we set a reduction factor of 2; this would give a comm size on c' of 2.
12278439623fSPatrick Sanan    Using the interlaced creation routine, the ranks in c with color = 0 will be rank 0 and 2.
122800fea0ebSDave May    We perform the scatter to the sub-communicator in the following way.
1229514db83aSDave May    [1] Given a vector x defined on communicator c
12306fc41876SBarry Smith 
1231514db83aSDave May .vb
1232514db83aSDave May    rank(c)  local values of x
1233514db83aSDave May    ------- ----------------------------------------
1234514db83aSDave May         0   [  0.0,  1.0,  2.0,  3.0,  4.0,  5.0 ]
1235514db83aSDave May         1   [  6.0,  7.0,  8.0,  9.0, 10.0, 11.0 ]
1236514db83aSDave May         2   [ 12.0, 13.0, 14.0, 15.0, 16.0, 17.0 ]
1237514db83aSDave May         3   [ 18.0, 19.0, 20.0, 21.0, 22.0, 23.0 ]
1238514db83aSDave May .ve
12396fc41876SBarry Smith 
1240514db83aSDave May    scatter into xtmp defined also on comm c, so that we have the following values
12416fc41876SBarry Smith 
1242514db83aSDave May .vb
1243514db83aSDave May    rank(c)  local values of xtmp
1244514db83aSDave May    ------- ----------------------------------------------------------------------------
1245514db83aSDave May         0   [  0.0,  1.0,  2.0,  3.0,  4.0,  5.0,  6.0,  7.0,  8.0,  9.0, 10.0, 11.0 ]
1246514db83aSDave May         1   [ ]
1247514db83aSDave May         2   [ 12.0, 13.0, 14.0, 15.0, 16.0, 17.0, 18.0, 19.0, 20.0, 21.0, 22.0, 23.0 ]
1248514db83aSDave May         3   [ ]
1249514db83aSDave May .ve
12506fc41876SBarry Smith 
12516fc41876SBarry Smith    The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') have no values
12526fc41876SBarry Smith 
1253514db83aSDave May    [2] Copy the values from ranks 0, 2 (indices with respect to comm c) into the vector xred which is defined on communicator c'.
12546fc41876SBarry Smith    Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0.
12556fc41876SBarry Smith 
1256514db83aSDave May .vb
1257514db83aSDave May    rank(c')  local values of xred
1258514db83aSDave May    -------- ----------------------------------------------------------------------------
1259514db83aSDave May          0   [  0.0,  1.0,  2.0,  3.0,  4.0,  5.0,  6.0,  7.0,  8.0,  9.0, 10.0, 11.0 ]
1260514db83aSDave May          1   [ 12.0, 13.0, 14.0, 15.0, 16.0, 17.0, 18.0, 19.0, 20.0, 21.0, 22.0, 23.0 ]
1261514db83aSDave May .ve
12621e07b27eSBarry Smith 
12631e07b27eSBarry Smith   Contributed by Dave May
12641e07b27eSBarry Smith 
1265bf00f589SPatrick Sanan   Reference:
1266bf00f589SPatrick Sanan   Dave A. May, Patrick Sanan, Karl Rupp, Matthew G. Knepley, and Barry F. Smith, "Extreme-Scale Multigrid Components within PETSc". 2016. In Proceedings of the Platform for Advanced Scientific Computing Conference (PASC '16). DOI: 10.1145/2929908.2929913
1267bf00f589SPatrick Sanan 
1268db781477SPatrick Sanan .seealso: `PCTelescopeGetKSP()`, `PCTelescopeGetDM()`, `PCTelescopeGetReductionFactor()`, `PCTelescopeSetReductionFactor()`, `PCTelescopeGetIgnoreDM()`, `PCTelescopeSetIgnoreDM()`, `PCREDUNDANT`
12691e07b27eSBarry Smith M*/
12709371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc) {
12711e07b27eSBarry Smith   struct _PC_Telescope *sred;
12721e07b27eSBarry Smith 
12731e07b27eSBarry Smith   PetscFunctionBegin;
12749566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc, &sred));
12752a22aa42SDave May   sred->psubcomm                   = NULL;
127648a10b22SPatrick Sanan   sred->subcommtype                = PETSC_SUBCOMM_INTERLACED;
12772a22aa42SDave May   sred->subcomm                    = MPI_COMM_NULL;
12781e07b27eSBarry Smith   sred->redfactor                  = 1;
12791e07b27eSBarry Smith   sred->ignore_dm                  = PETSC_FALSE;
12807c5279cbSDave May   sred->ignore_kspcomputeoperators = PETSC_FALSE;
12818d9f7141SDave May   sred->use_coarse_dm              = PETSC_FALSE;
12821e07b27eSBarry Smith   pc->data                         = (void *)sred;
12831e07b27eSBarry Smith 
12841e07b27eSBarry Smith   pc->ops->apply           = PCApply_Telescope;
12851e07b27eSBarry Smith   pc->ops->applytranspose  = NULL;
1286f650675bSDave May   pc->ops->applyrichardson = PCApplyRichardson_Telescope;
12871e07b27eSBarry Smith   pc->ops->setup           = PCSetUp_Telescope;
12881e07b27eSBarry Smith   pc->ops->destroy         = PCDestroy_Telescope;
12891e07b27eSBarry Smith   pc->ops->reset           = PCReset_Telescope;
12901e07b27eSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Telescope;
12911e07b27eSBarry Smith   pc->ops->view            = PCView_Telescope;
12921e07b27eSBarry Smith 
12931e07b27eSBarry Smith   sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
12941e07b27eSBarry Smith   sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
12951e07b27eSBarry Smith   sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
12961e07b27eSBarry Smith   sred->pctelescope_reset_type              = NULL;
12971e07b27eSBarry Smith 
12989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetKSP_C", PCTelescopeGetKSP_Telescope));
12999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetSubcommType_C", PCTelescopeGetSubcommType_Telescope));
13009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetSubcommType_C", PCTelescopeSetSubcommType_Telescope));
13019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetReductionFactor_C", PCTelescopeGetReductionFactor_Telescope));
13029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetReductionFactor_C", PCTelescopeSetReductionFactor_Telescope));
13039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetIgnoreDM_C", PCTelescopeGetIgnoreDM_Telescope));
13049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetIgnoreDM_C", PCTelescopeSetIgnoreDM_Telescope));
13059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetIgnoreKSPComputeOperators_C", PCTelescopeGetIgnoreKSPComputeOperators_Telescope));
13069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetIgnoreKSPComputeOperators_C", PCTelescopeSetIgnoreKSPComputeOperators_Telescope));
13079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetDM_C", PCTelescopeGetDM_Telescope));
13089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetUseCoarseDM_C", PCTelescopeGetUseCoarseDM_Telescope));
13099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetUseCoarseDM_C", PCTelescopeSetUseCoarseDM_Telescope));
13101e07b27eSBarry Smith   PetscFunctionReturn(0);
13111e07b27eSBarry Smith }
1312