xref: /petsc/src/ksp/pc/impls/telescope/telescope.c (revision b5f0bcd6e9e8ed97648738542f5163d94f7b1782)
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 */
PCTelescopeTestValidSubcomm(MPI_Comm comm_f,MPI_Comm comm_c,PetscBool * isvalid)5166976f2fSJacob Faibussowitsch static PetscErrorCode PCTelescopeTestValidSubcomm(MPI_Comm comm_f, MPI_Comm comm_c, PetscBool *isvalid)
52d71ae5a4SJacob Faibussowitsch {
5357f12427SDave May   PetscInt     valid = 1;
548d9f7141SDave May   MPI_Group    group_f, group_c;
558d9f7141SDave May   PetscMPIInt  count, k, size_f = 0, size_c = 0, size_c_sum = 0;
565bd1e576SStefano Zampini   PetscMPIInt *ranks_f, *ranks_c;
578d9f7141SDave May 
5857f12427SDave May   PetscFunctionBegin;
5908401ef6SPierre Jolivet   PetscCheck(comm_f != MPI_COMM_NULL, PETSC_COMM_SELF, PETSC_ERR_SUP, "comm_f cannot be MPI_COMM_NULL");
608d9f7141SDave May 
619566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_group(comm_f, &group_f));
6248a46eb9SPierre Jolivet   if (comm_c != MPI_COMM_NULL) PetscCallMPI(MPI_Comm_group(comm_c, &group_c));
638d9f7141SDave May 
649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm_f, &size_f));
6548a46eb9SPierre Jolivet   if (comm_c != MPI_COMM_NULL) PetscCallMPI(MPI_Comm_size(comm_c, &size_c));
668d9f7141SDave May 
678d9f7141SDave May   /* check not all comm_c's are NULL */
688d9f7141SDave May   size_c_sum = size_c;
696a210b70SBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &size_c_sum, 1, MPI_INT, MPI_SUM, comm_f));
705bd1e576SStefano Zampini   if (size_c_sum == 0) valid = 0;
718d9f7141SDave May 
728d9f7141SDave May   /* check we can map at least 1 rank in comm_c to comm_f */
739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size_f, &ranks_f));
749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size_c, &ranks_c));
755bd1e576SStefano Zampini   for (k = 0; k < size_f; k++) ranks_f[k] = MPI_UNDEFINED;
765bd1e576SStefano Zampini   for (k = 0; k < size_c; k++) ranks_c[k] = k;
778d9f7141SDave May 
7857f12427SDave May   /*
7957f12427SDave May    MPI_Group_translate_ranks() returns a non-zero exit code if any rank cannot be translated.
8057f12427SDave May    I do not want the code to terminate immediately if this occurs, rather I want to throw
8157f12427SDave May    the error later (during PCSetUp_Telescope()) via SETERRQ() with a message indicating
8257f12427SDave May    that comm_c is not a valid sub-communicator.
839566063dSJacob Faibussowitsch    Hence I purposefully do not call PetscCall() after MPI_Group_translate_ranks().
8457f12427SDave May   */
858d9f7141SDave May   count = 0;
868d9f7141SDave May   if (comm_c != MPI_COMM_NULL) {
8766b79024SDave May     (void)MPI_Group_translate_ranks(group_c, size_c, ranks_c, group_f, ranks_f);
888d9f7141SDave May     for (k = 0; k < size_f; k++) {
89ad540459SPierre Jolivet       if (ranks_f[k] == MPI_UNDEFINED) count++;
908d9f7141SDave May     }
918d9f7141SDave May   }
925bd1e576SStefano Zampini   if (count == size_f) valid = 0;
938d9f7141SDave May 
94462c564dSBarry Smith   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &valid, 1, MPIU_INT, MPI_MIN, comm_f));
955bd1e576SStefano Zampini   if (valid == 1) *isvalid = PETSC_TRUE;
965bd1e576SStefano Zampini   else *isvalid = PETSC_FALSE;
978d9f7141SDave May 
989566063dSJacob Faibussowitsch   PetscCall(PetscFree(ranks_f));
999566063dSJacob Faibussowitsch   PetscCall(PetscFree(ranks_c));
1009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Group_free(&group_f));
10148a46eb9SPierre Jolivet   if (comm_c != MPI_COMM_NULL) PetscCallMPI(MPI_Group_free(&group_c));
1023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1038d9f7141SDave May }
1048d9f7141SDave May 
private_PCTelescopeGetSubDM(PC_Telescope sred)10566976f2fSJacob Faibussowitsch static DM private_PCTelescopeGetSubDM(PC_Telescope sred)
106d71ae5a4SJacob Faibussowitsch {
107c6a0d831SBarry Smith   DM subdm = NULL;
1081e07b27eSBarry Smith 
1099371c9d4SSatish Balay   if (!PCTelescope_isActiveRank(sred)) {
1109371c9d4SSatish Balay     subdm = NULL;
1119371c9d4SSatish Balay   } else {
1121e07b27eSBarry Smith     switch (sred->sr_type) {
113d71ae5a4SJacob Faibussowitsch     case TELESCOPE_DEFAULT:
114d71ae5a4SJacob Faibussowitsch       subdm = NULL;
115d71ae5a4SJacob Faibussowitsch       break;
116d71ae5a4SJacob Faibussowitsch     case TELESCOPE_DMDA:
117d71ae5a4SJacob Faibussowitsch       subdm = ((PC_Telescope_DMDACtx *)sred->dm_ctx)->dmrepart;
118d71ae5a4SJacob Faibussowitsch       break;
119d71ae5a4SJacob Faibussowitsch     case TELESCOPE_DMPLEX:
120d71ae5a4SJacob Faibussowitsch       subdm = NULL;
121d71ae5a4SJacob Faibussowitsch       break;
1229371c9d4SSatish Balay     case TELESCOPE_COARSEDM:
1233ba16761SJacob Faibussowitsch       if (sred->ksp) PetscCallAbort(PETSC_COMM_SELF, KSPGetDM(sred->ksp, &subdm));
1248d9f7141SDave May       break;
1251e07b27eSBarry Smith     }
1261e07b27eSBarry Smith   }
1273ba16761SJacob Faibussowitsch   return subdm;
1281e07b27eSBarry Smith }
1291e07b27eSBarry Smith 
PCTelescopeSetUp_default(PC pc,PC_Telescope sred)13066976f2fSJacob Faibussowitsch static PetscErrorCode PCTelescopeSetUp_default(PC pc, PC_Telescope sred)
131d71ae5a4SJacob Faibussowitsch {
1321e07b27eSBarry Smith   PetscInt   m, M, bs, st, ed;
1331e07b27eSBarry Smith   Vec        x, xred, yred, xtmp;
1341e07b27eSBarry Smith   Mat        B;
1351e07b27eSBarry Smith   MPI_Comm   comm, subcomm;
1361e07b27eSBarry Smith   VecScatter scatter;
1371e07b27eSBarry Smith   IS         isin;
13854cd43dcSJunchao Zhang   VecType    vectype;
1391e07b27eSBarry Smith 
1401e07b27eSBarry Smith   PetscFunctionBegin;
1419566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "PCTelescope: setup (default)\n"));
1421e07b27eSBarry Smith   comm    = PetscSubcommParent(sred->psubcomm);
1431e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
1441e07b27eSBarry Smith 
1459566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &B));
1469566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, &M, NULL));
1479566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(B, &bs));
1489566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(B, &x, NULL));
1497addb90fSBarry Smith   PetscCall(MatGetVecType(B, &vectype)); /* Use the vectype of the matrix used to construct the preconditioner by default */
1501e07b27eSBarry Smith 
1511e07b27eSBarry Smith   xred = NULL;
1523ac26c5eSBarry Smith   m    = 0;
15357f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
1549566063dSJacob Faibussowitsch     PetscCall(VecCreate(subcomm, &xred));
1559566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(xred, PETSC_DECIDE, M));
1569566063dSJacob Faibussowitsch     PetscCall(VecSetBlockSize(xred, bs));
1577addb90fSBarry Smith     PetscCall(VecSetType(xred, vectype));
1589566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(xred));
1599566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(xred, &m));
1601e07b27eSBarry Smith   }
1611e07b27eSBarry Smith 
1621e07b27eSBarry Smith   yred = NULL;
16348a46eb9SPierre Jolivet   if (PCTelescope_isActiveRank(sred)) PetscCall(VecDuplicate(xred, &yred));
1641e07b27eSBarry Smith 
1659566063dSJacob Faibussowitsch   PetscCall(VecCreate(comm, &xtmp));
1669566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(xtmp, m, PETSC_DECIDE));
1679566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(xtmp, bs));
1689566063dSJacob Faibussowitsch   PetscCall(VecSetType(xtmp, vectype));
1691e07b27eSBarry Smith 
17057f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
1719566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(xred, &st, &ed));
17257508eceSPierre Jolivet     PetscCall(ISCreateStride(comm, ed - st, st, 1, &isin));
1731e07b27eSBarry Smith   } else {
1749566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(x, &st, &ed));
1759566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, 0, st, 1, &isin));
1761e07b27eSBarry Smith   }
1779566063dSJacob Faibussowitsch   PetscCall(ISSetBlockSize(isin, bs));
1781e07b27eSBarry Smith 
1799566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(x, isin, xtmp, NULL, &scatter));
1801e07b27eSBarry Smith 
1811e07b27eSBarry Smith   sred->isin    = isin;
1821e07b27eSBarry Smith   sred->scatter = scatter;
1831e07b27eSBarry Smith   sred->xred    = xred;
1841e07b27eSBarry Smith   sred->yred    = yred;
1851e07b27eSBarry Smith   sred->xtmp    = xtmp;
1869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1881e07b27eSBarry Smith }
1891e07b27eSBarry Smith 
PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat * A)19066976f2fSJacob Faibussowitsch static PetscErrorCode PCTelescopeMatCreate_default(PC pc, PC_Telescope sred, MatReuse reuse, Mat *A)
191d71ae5a4SJacob Faibussowitsch {
1921e07b27eSBarry Smith   MPI_Comm comm, subcomm;
1931e07b27eSBarry Smith   Mat      Bred, B;
1945624f943SPierre Jolivet   PetscInt nr, nc, bs;
1951e07b27eSBarry Smith   IS       isrow, iscol;
1961e07b27eSBarry Smith   Mat      Blocal, *_Blocal;
1971e07b27eSBarry Smith 
1981e07b27eSBarry Smith   PetscFunctionBegin;
1999566063dSJacob Faibussowitsch   PetscCall(PetscInfo(pc, "PCTelescope: updating the redundant preconditioned operator (default)\n"));
2009566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
2011e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
2029566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &B));
2039566063dSJacob Faibussowitsch   PetscCall(MatGetSize(B, &nr, &nc));
2041e07b27eSBarry Smith   isrow = sred->isin;
2059566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, nc, 0, 1, &iscol));
2069566063dSJacob Faibussowitsch   PetscCall(ISSetIdentity(iscol));
2079566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSizes(B, NULL, &bs));
2089566063dSJacob Faibussowitsch   PetscCall(ISSetBlockSize(iscol, bs));
2099566063dSJacob Faibussowitsch   PetscCall(MatSetOption(B, MAT_SUBMAT_SINGLEIS, PETSC_TRUE));
2109566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, 1, &isrow, &iscol, MAT_INITIAL_MATRIX, &_Blocal));
2111e07b27eSBarry Smith   Blocal = *_Blocal;
2129566063dSJacob Faibussowitsch   PetscCall(PetscFree(_Blocal));
2131e07b27eSBarry Smith   Bred = NULL;
21457f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
2151e07b27eSBarry Smith     PetscInt mm;
2161e07b27eSBarry Smith 
217ad540459SPierre Jolivet     if (reuse != MAT_INITIAL_MATRIX) Bred = *A;
2181e07b27eSBarry Smith 
2199566063dSJacob Faibussowitsch     PetscCall(MatGetSize(Blocal, &mm, NULL));
2209566063dSJacob Faibussowitsch     PetscCall(MatCreateMPIMatConcatenateSeqMat(subcomm, Blocal, mm, reuse, &Bred));
2219a26306aSPierre Jolivet     PetscCall(MatPropagateSymmetryOptions(B, Bred));
2221e07b27eSBarry Smith   }
2231e07b27eSBarry Smith   *A = Bred;
2249566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscol));
2259566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Blocal));
2263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2271e07b27eSBarry Smith }
2281e07b27eSBarry Smith 
PCTelescopeSubNullSpaceCreate_Telescope(PC pc,PC_Telescope sred,MatNullSpace nullspace,MatNullSpace * sub_nullspace)229d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeSubNullSpaceCreate_Telescope(PC pc, PC_Telescope sred, MatNullSpace nullspace, MatNullSpace *sub_nullspace)
230d71ae5a4SJacob Faibussowitsch {
2311e07b27eSBarry Smith   PetscBool  has_const;
2321e07b27eSBarry Smith   const Vec *vecs;
233c41e779fSDave May   Vec       *sub_vecs = NULL;
234392968a1SPatrick Sanan   PetscInt   i, k, n = 0;
2351e07b27eSBarry Smith   MPI_Comm   subcomm;
2361e07b27eSBarry Smith 
2371e07b27eSBarry Smith   PetscFunctionBegin;
2381e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
2399566063dSJacob Faibussowitsch   PetscCall(MatNullSpaceGetVecs(nullspace, &has_const, &n, &vecs));
2401e07b27eSBarry Smith 
24157f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
24248a46eb9SPierre Jolivet     if (n) PetscCall(VecDuplicateVecs(sred->xred, n, &sub_vecs));
2431e07b27eSBarry Smith   }
2441e07b27eSBarry Smith 
2451e07b27eSBarry Smith   /* copy entries */
2461e07b27eSBarry Smith   for (k = 0; k < n; k++) {
2471e07b27eSBarry Smith     const PetscScalar *x_array;
2481e07b27eSBarry Smith     PetscScalar       *LA_sub_vec;
24913c30530SDave May     PetscInt           st, ed;
2501e07b27eSBarry Smith 
2511e07b27eSBarry Smith     /* pull in vector x->xtmp */
2529566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(sred->scatter, vecs[k], sred->xtmp, INSERT_VALUES, SCATTER_FORWARD));
2539566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(sred->scatter, vecs[k], sred->xtmp, INSERT_VALUES, SCATTER_FORWARD));
25447856c66SBarry Smith     if (sub_vecs) {
255a04a6428SPatrick Sanan       /* copy vector entries into xred */
2569566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(sred->xtmp, &x_array));
257ea2b237eSDave May       if (sub_vecs[k]) {
2589566063dSJacob Faibussowitsch         PetscCall(VecGetOwnershipRange(sub_vecs[k], &st, &ed));
2599566063dSJacob Faibussowitsch         PetscCall(VecGetArray(sub_vecs[k], &LA_sub_vec));
260ad540459SPierre Jolivet         for (i = 0; i < ed - st; i++) LA_sub_vec[i] = x_array[i];
2619566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(sub_vecs[k], &LA_sub_vec));
2621e07b27eSBarry Smith       }
2639566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(sred->xtmp, &x_array));
2641e07b27eSBarry Smith     }
26547856c66SBarry Smith   }
2661e07b27eSBarry Smith 
26757f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
268d8b9d5b7SPatrick Sanan     /* create new (near) nullspace for redundant object */
2699566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceCreate(subcomm, has_const, n, sub_vecs, sub_nullspace));
2709566063dSJacob Faibussowitsch     PetscCall(VecDestroyVecs(n, &sub_vecs));
27128b400f6SJacob Faibussowitsch     PetscCheck(!nullspace->remove, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope");
27228b400f6SJacob Faibussowitsch     PetscCheck(!nullspace->rmctx, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Propagation of custom remove callback context not supported when propagating (near) nullspaces with PCTelescope");
273d8b9d5b7SPatrick Sanan   }
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
275392968a1SPatrick Sanan }
276392968a1SPatrick Sanan 
PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat)277d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc, PC_Telescope sred, Mat sub_mat)
278d71ae5a4SJacob Faibussowitsch {
279392968a1SPatrick Sanan   Mat B;
280392968a1SPatrick Sanan 
281392968a1SPatrick Sanan   PetscFunctionBegin;
2829566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc, NULL, &B));
283392968a1SPatrick Sanan   /* Propagate the nullspace if it exists */
284392968a1SPatrick Sanan   {
285392968a1SPatrick Sanan     MatNullSpace nullspace, sub_nullspace;
2869566063dSJacob Faibussowitsch     PetscCall(MatGetNullSpace(B, &nullspace));
287392968a1SPatrick Sanan     if (nullspace) {
2889566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "PCTelescope: generating nullspace (default)\n"));
2899566063dSJacob Faibussowitsch       PetscCall(PCTelescopeSubNullSpaceCreate_Telescope(pc, sred, nullspace, &sub_nullspace));
29057f12427SDave May       if (PCTelescope_isActiveRank(sred)) {
2919566063dSJacob Faibussowitsch         PetscCall(MatSetNullSpace(sub_mat, sub_nullspace));
2929566063dSJacob Faibussowitsch         PetscCall(MatNullSpaceDestroy(&sub_nullspace));
2931e07b27eSBarry Smith       }
294392968a1SPatrick Sanan     }
295392968a1SPatrick Sanan   }
296392968a1SPatrick Sanan   /* Propagate the near nullspace if it exists */
297392968a1SPatrick Sanan   {
298392968a1SPatrick Sanan     MatNullSpace nearnullspace, sub_nearnullspace;
2999566063dSJacob Faibussowitsch     PetscCall(MatGetNearNullSpace(B, &nearnullspace));
300392968a1SPatrick Sanan     if (nearnullspace) {
3019566063dSJacob Faibussowitsch       PetscCall(PetscInfo(pc, "PCTelescope: generating near nullspace (default)\n"));
3029566063dSJacob Faibussowitsch       PetscCall(PCTelescopeSubNullSpaceCreate_Telescope(pc, sred, nearnullspace, &sub_nearnullspace));
30357f12427SDave May       if (PCTelescope_isActiveRank(sred)) {
3049566063dSJacob Faibussowitsch         PetscCall(MatSetNearNullSpace(sub_mat, sub_nearnullspace));
3059566063dSJacob Faibussowitsch         PetscCall(MatNullSpaceDestroy(&sub_nearnullspace));
306392968a1SPatrick Sanan       }
307392968a1SPatrick Sanan     }
308392968a1SPatrick Sanan   }
3093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3101e07b27eSBarry Smith }
3111e07b27eSBarry Smith 
PCView_Telescope(PC pc,PetscViewer viewer)312d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Telescope(PC pc, PetscViewer viewer)
313d71ae5a4SJacob Faibussowitsch {
3141e07b27eSBarry Smith   PC_Telescope sred = (PC_Telescope)pc->data;
3159f196a02SMartin Diehl   PetscBool    isascii, isstring;
3161e07b27eSBarry Smith   PetscViewer  subviewer;
3171e07b27eSBarry Smith 
3181e07b27eSBarry Smith   PetscFunctionBegin;
3199f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
3209566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
3219f196a02SMartin Diehl   if (isascii) {
3228d9f7141SDave May     {
3231e07b27eSBarry Smith       MPI_Comm    comm, subcomm;
3241e07b27eSBarry Smith       PetscMPIInt comm_size, subcomm_size;
3258d9f7141SDave May       DM          dm = NULL, subdm = NULL;
3261e07b27eSBarry Smith 
3279566063dSJacob Faibussowitsch       PetscCall(PCGetDM(pc, &dm));
3281e07b27eSBarry Smith       subdm = private_PCTelescopeGetSubDM(sred);
3298d9f7141SDave May 
3308d9f7141SDave May       if (sred->psubcomm) {
3311e07b27eSBarry Smith         comm    = PetscSubcommParent(sred->psubcomm);
3321e07b27eSBarry Smith         subcomm = PetscSubcommChild(sred->psubcomm);
3339566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(comm, &comm_size));
3349566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(subcomm, &subcomm_size));
3351e07b27eSBarry Smith 
3369566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
337f0b74427SPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "PETSc subcomm: parent comm size reduction factor = %" PetscInt_FMT "\n", sred->redfactor));
338f0b74427SPierre Jolivet         PetscCall(PetscViewerASCIIPrintf(viewer, "PETSc subcomm: parent_size = %d , subcomm_size = %d\n", comm_size, subcomm_size));
33948a10b22SPatrick Sanan         switch (sred->subcommtype) {
340d71ae5a4SJacob Faibussowitsch         case PETSC_SUBCOMM_INTERLACED:
341f0b74427SPierre Jolivet           PetscCall(PetscViewerASCIIPrintf(viewer, "PETSc subcomm: type = %s\n", PetscSubcommTypes[sred->subcommtype]));
342d71ae5a4SJacob Faibussowitsch           break;
343d71ae5a4SJacob Faibussowitsch         case PETSC_SUBCOMM_CONTIGUOUS:
344f0b74427SPierre Jolivet           PetscCall(PetscViewerASCIIPrintf(viewer, "PETSc subcomm type = %s\n", PetscSubcommTypes[sred->subcommtype]));
345d71ae5a4SJacob Faibussowitsch           break;
346d71ae5a4SJacob Faibussowitsch         default:
347d71ae5a4SJacob Faibussowitsch           SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "General subcomm type not supported by PCTelescope");
34848a10b22SPatrick Sanan         }
3499566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
3508d9f7141SDave May       } else {
3519566063dSJacob Faibussowitsch         PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
3528d9f7141SDave May         subcomm = sred->subcomm;
353ad540459SPierre Jolivet         if (!PCTelescope_isActiveRank(sred)) subcomm = PETSC_COMM_SELF;
3548d9f7141SDave May 
3559566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
3569566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "subcomm: using user provided sub-communicator\n"));
3579566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
3588d9f7141SDave May       }
3598d9f7141SDave May 
3609566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, subcomm, &subviewer));
36157f12427SDave May       if (PCTelescope_isActiveRank(sred)) {
3629566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(subviewer));
3631e07b27eSBarry Smith 
36448a46eb9SPierre Jolivet         if (dm && sred->ignore_dm) PetscCall(PetscViewerASCIIPrintf(subviewer, "ignoring DM\n"));
36548a46eb9SPierre Jolivet         if (sred->ignore_kspcomputeoperators) PetscCall(PetscViewerASCIIPrintf(subviewer, "ignoring KSPComputeOperators\n"));
3661e07b27eSBarry Smith         switch (sred->sr_type) {
367d71ae5a4SJacob Faibussowitsch         case TELESCOPE_DEFAULT:
368d71ae5a4SJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "setup type: default\n"));
369d71ae5a4SJacob Faibussowitsch           break;
3701e07b27eSBarry Smith         case TELESCOPE_DMDA:
3719566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "setup type: DMDA auto-repartitioning\n"));
3729566063dSJacob Faibussowitsch           PetscCall(DMView_DA_Short(subdm, subviewer));
3731e07b27eSBarry Smith           break;
374d71ae5a4SJacob Faibussowitsch         case TELESCOPE_DMPLEX:
375d71ae5a4SJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "setup type: DMPLEX auto-repartitioning\n"));
376d71ae5a4SJacob Faibussowitsch           break;
377d71ae5a4SJacob Faibussowitsch         case TELESCOPE_COARSEDM:
378d71ae5a4SJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "setup type: coarse DM\n"));
379d71ae5a4SJacob Faibussowitsch           break;
3808d9f7141SDave May         }
3818d9f7141SDave May 
3828d9f7141SDave May         if (dm) {
3838d9f7141SDave May           PetscObject obj = (PetscObject)dm;
3849566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "Parent DM object:"));
3853ba16761SJacob Faibussowitsch           PetscCall(PetscViewerASCIIUseTabs(subviewer, PETSC_FALSE));
3863ba16761SJacob Faibussowitsch           if (obj->type_name) PetscCall(PetscViewerASCIIPrintf(subviewer, " type = %s;", obj->type_name));
3873ba16761SJacob Faibussowitsch           if (obj->name) PetscCall(PetscViewerASCIIPrintf(subviewer, " name = %s;", obj->name));
3883ba16761SJacob Faibussowitsch           if (obj->prefix) PetscCall(PetscViewerASCIIPrintf(subviewer, " prefix = %s", obj->prefix));
3899566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "\n"));
3903ba16761SJacob Faibussowitsch           PetscCall(PetscViewerASCIIUseTabs(subviewer, PETSC_TRUE));
3918d9f7141SDave May         } else {
3929566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "Parent DM object: NULL\n"));
3938d9f7141SDave May         }
3948d9f7141SDave May         if (subdm) {
3958d9f7141SDave May           PetscObject obj = (PetscObject)subdm;
3969566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "Sub DM object:"));
3973ba16761SJacob Faibussowitsch           PetscCall(PetscViewerASCIIUseTabs(subviewer, PETSC_FALSE));
3983ba16761SJacob Faibussowitsch           if (obj->type_name) PetscCall(PetscViewerASCIIPrintf(subviewer, " type = %s;", obj->type_name));
3993ba16761SJacob Faibussowitsch           if (obj->name) PetscCall(PetscViewerASCIIPrintf(subviewer, " name = %s;", obj->name));
4003ba16761SJacob Faibussowitsch           if (obj->prefix) PetscCall(PetscViewerASCIIPrintf(subviewer, " prefix = %s", obj->prefix));
4019566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "\n"));
4023ba16761SJacob Faibussowitsch           PetscCall(PetscViewerASCIIUseTabs(subviewer, PETSC_TRUE));
4038d9f7141SDave May         } else {
4049566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPrintf(subviewer, "Sub DM object: NULL\n"));
4051e07b27eSBarry Smith         }
4061e07b27eSBarry Smith 
4079566063dSJacob Faibussowitsch         PetscCall(KSPView(sred->ksp, subviewer));
4089566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(subviewer));
4091e07b27eSBarry Smith       }
4109566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, subcomm, &subviewer));
4111e07b27eSBarry Smith     }
4121e07b27eSBarry Smith   }
4133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4141e07b27eSBarry Smith }
4151e07b27eSBarry Smith 
PCSetUp_Telescope(PC pc)416d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Telescope(PC pc)
417d71ae5a4SJacob Faibussowitsch {
4181e07b27eSBarry Smith   PC_Telescope    sred = (PC_Telescope)pc->data;
419bd49479cSSatish Balay   MPI_Comm        comm, subcomm = 0;
4201e07b27eSBarry Smith   PCTelescopeType sr_type;
4211e07b27eSBarry Smith 
4221e07b27eSBarry Smith   PetscFunctionBegin;
4239566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
4241e07b27eSBarry Smith 
4251e07b27eSBarry Smith   /* Determine type of setup/update */
4261e07b27eSBarry Smith   if (!pc->setupcalled) {
4271e07b27eSBarry Smith     PetscBool has_dm, same;
4281e07b27eSBarry Smith     DM        dm;
4291e07b27eSBarry Smith 
4301e07b27eSBarry Smith     sr_type = TELESCOPE_DEFAULT;
4311e07b27eSBarry Smith     has_dm  = PETSC_FALSE;
4329566063dSJacob Faibussowitsch     PetscCall(PCGetDM(pc, &dm));
433ad540459SPierre Jolivet     if (dm) has_dm = PETSC_TRUE;
4341e07b27eSBarry Smith     if (has_dm) {
4351e07b27eSBarry Smith       /* check for dmda */
4369566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMDA, &same));
4371e07b27eSBarry Smith       if (same) {
4389566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "PCTelescope: found DMDA\n"));
4391e07b27eSBarry Smith         sr_type = TELESCOPE_DMDA;
4401e07b27eSBarry Smith       }
4411e07b27eSBarry Smith       /* check for dmplex */
4429566063dSJacob Faibussowitsch       PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMPLEX, &same));
4431e07b27eSBarry Smith       if (same) {
4449566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "PCTelescope: found DMPLEX\n"));
4451e07b27eSBarry Smith         sr_type = TELESCOPE_DMPLEX;
4461e07b27eSBarry Smith       }
4478d9f7141SDave May 
4488d9f7141SDave May       if (sred->use_coarse_dm) {
4499566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "PCTelescope: using coarse DM\n"));
4508d9f7141SDave May         sr_type = TELESCOPE_COARSEDM;
4511e07b27eSBarry Smith       }
4521e07b27eSBarry Smith 
4531e07b27eSBarry Smith       if (sred->ignore_dm) {
4549566063dSJacob Faibussowitsch         PetscCall(PetscInfo(pc, "PCTelescope: ignoring DM\n"));
4551e07b27eSBarry Smith         sr_type = TELESCOPE_DEFAULT;
4561e07b27eSBarry Smith       }
4578d9f7141SDave May     }
4581e07b27eSBarry Smith     sred->sr_type = sr_type;
4591e07b27eSBarry Smith   } else {
4601e07b27eSBarry Smith     sr_type = sred->sr_type;
4611e07b27eSBarry Smith   }
4621e07b27eSBarry Smith 
463d8b9d5b7SPatrick Sanan   /* set function pointers for repartition setup, matrix creation/update, matrix (near) nullspace, and reset functionality */
4641e07b27eSBarry Smith   switch (sr_type) {
4651e07b27eSBarry Smith   case TELESCOPE_DEFAULT:
4661e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
4671e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
4681e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
4691e07b27eSBarry Smith     sred->pctelescope_reset_type              = NULL;
4701e07b27eSBarry Smith     break;
4711e07b27eSBarry Smith   case TELESCOPE_DMDA:
4721e07b27eSBarry Smith     pc->ops->apply                            = PCApply_Telescope_dmda;
473f650675bSDave May     pc->ops->applyrichardson                  = PCApplyRichardson_Telescope_dmda;
4741e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_dmda;
4751e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_dmda;
4761e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda;
4771e07b27eSBarry Smith     sred->pctelescope_reset_type              = PCReset_Telescope_dmda;
4781e07b27eSBarry Smith     break;
479d71ae5a4SJacob Faibussowitsch   case TELESCOPE_DMPLEX:
480d71ae5a4SJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_SUP, "Support for DMPLEX is currently not available");
4818d9f7141SDave May   case TELESCOPE_COARSEDM:
4828d9f7141SDave May     pc->ops->apply                            = PCApply_Telescope_CoarseDM;
4838d9f7141SDave May     pc->ops->applyrichardson                  = PCApplyRichardson_Telescope_CoarseDM;
4848d9f7141SDave May     sred->pctelescope_setup_type              = PCTelescopeSetUp_CoarseDM;
4858d9f7141SDave May     sred->pctelescope_matcreate_type          = NULL;
4868d9f7141SDave May     sred->pctelescope_matnullspacecreate_type = NULL; /* PCTelescopeMatNullSpaceCreate_CoarseDM; */
4878d9f7141SDave May     sred->pctelescope_reset_type              = PCReset_Telescope_CoarseDM;
4881e07b27eSBarry Smith     break;
489d71ae5a4SJacob Faibussowitsch   default:
490d71ae5a4SJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_SUP, "Support only provided for: repartitioning an operator; repartitioning a DMDA; or using a coarse DM");
4918d9f7141SDave May   }
4928d9f7141SDave May 
4938d9f7141SDave May   /* subcomm definition */
4948d9f7141SDave May   if (!pc->setupcalled) {
4958d9f7141SDave May     if ((sr_type == TELESCOPE_DEFAULT) || (sr_type == TELESCOPE_DMDA)) {
4968d9f7141SDave May       if (!sred->psubcomm) {
4979566063dSJacob Faibussowitsch         PetscCall(PetscSubcommCreate(comm, &sred->psubcomm));
4989566063dSJacob Faibussowitsch         PetscCall(PetscSubcommSetNumber(sred->psubcomm, sred->redfactor));
4999566063dSJacob Faibussowitsch         PetscCall(PetscSubcommSetType(sred->psubcomm, sred->subcommtype));
5008d9f7141SDave May         sred->subcomm = PetscSubcommChild(sred->psubcomm);
5018d9f7141SDave May       }
5028d9f7141SDave May     } else { /* query PC for DM, check communicators */
5038d9f7141SDave May       DM          dm, dm_coarse_partition          = NULL;
5048d9f7141SDave May       MPI_Comm    comm_fine, comm_coarse_partition = MPI_COMM_NULL;
5058d9f7141SDave May       PetscMPIInt csize_fine = 0, csize_coarse_partition = 0, cs[2], csg[2], cnt = 0;
50666976f2fSJacob Faibussowitsch       PetscBool   isvalidsubcomm = PETSC_TRUE;
5078d9f7141SDave May 
5089566063dSJacob Faibussowitsch       PetscCall(PCGetDM(pc, &dm));
5098d9f7141SDave May       comm_fine = PetscObjectComm((PetscObject)dm);
5109566063dSJacob Faibussowitsch       PetscCall(DMGetCoarseDM(dm, &dm_coarse_partition));
511ad540459SPierre Jolivet       if (dm_coarse_partition) cnt = 1;
5126a210b70SBarry Smith       PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &cnt, 1, MPI_INT, MPI_SUM, comm_fine));
51308401ef6SPierre Jolivet       PetscCheck(cnt != 0, comm_fine, PETSC_ERR_SUP, "Zero instances of a coarse DM were found");
5148d9f7141SDave May 
5159566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_size(comm_fine, &csize_fine));
5168d9f7141SDave May       if (dm_coarse_partition) {
5178d9f7141SDave May         comm_coarse_partition = PetscObjectComm((PetscObject)dm_coarse_partition);
5189566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(comm_coarse_partition, &csize_coarse_partition));
5198d9f7141SDave May       }
5208d9f7141SDave May 
5218d9f7141SDave May       cs[0] = csize_fine;
5228d9f7141SDave May       cs[1] = csize_coarse_partition;
5236a210b70SBarry Smith       PetscCallMPI(MPIU_Allreduce(cs, csg, 2, MPI_INT, MPI_MAX, comm_fine));
52408401ef6SPierre 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");
5258d9f7141SDave May 
5269566063dSJacob Faibussowitsch       PetscCall(PCTelescopeTestValidSubcomm(comm_fine, comm_coarse_partition, &isvalidsubcomm));
52728b400f6SJacob Faibussowitsch       PetscCheck(isvalidsubcomm, comm_fine, PETSC_ERR_SUP, "Coarse DM communicator is not a sub-communicator of parentDM->comm");
5288d9f7141SDave May       sred->subcomm = comm_coarse_partition;
5298d9f7141SDave May     }
5308d9f7141SDave May   }
5318d9f7141SDave May   subcomm = sred->subcomm;
5328d9f7141SDave May 
5338d9f7141SDave May   /* internal KSP */
5348d9f7141SDave May   if (!pc->setupcalled) {
5358d9f7141SDave May     const char *prefix;
5368d9f7141SDave May 
53757f12427SDave May     if (PCTelescope_isActiveRank(sred)) {
5389566063dSJacob Faibussowitsch       PetscCall(KSPCreate(subcomm, &sred->ksp));
5393821be0aSBarry Smith       PetscCall(KSPSetNestLevel(sred->ksp, pc->kspnestlevel));
5409566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(sred->ksp, pc->erroriffailure));
5419566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)sred->ksp, (PetscObject)pc, 1));
5429566063dSJacob Faibussowitsch       PetscCall(KSPSetType(sred->ksp, KSPPREONLY));
5439566063dSJacob Faibussowitsch       PetscCall(PCGetOptionsPrefix(pc, &prefix));
5449566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(sred->ksp, prefix));
5459566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(sred->ksp, "telescope_"));
5468d9f7141SDave May     }
5471e07b27eSBarry Smith   }
5481e07b27eSBarry Smith 
5491e07b27eSBarry Smith   /* setup */
55048a46eb9SPierre Jolivet   if (!pc->setupcalled && sred->pctelescope_setup_type) PetscCall(sred->pctelescope_setup_type(pc, sred));
5511e07b27eSBarry Smith   /* update */
5521e07b27eSBarry Smith   if (!pc->setupcalled) {
55348a46eb9SPierre Jolivet     if (sred->pctelescope_matcreate_type) PetscCall(sred->pctelescope_matcreate_type(pc, sred, MAT_INITIAL_MATRIX, &sred->Bred));
5541baa6e33SBarry Smith     if (sred->pctelescope_matnullspacecreate_type) PetscCall(sred->pctelescope_matnullspacecreate_type(pc, sred, sred->Bred));
5551e07b27eSBarry Smith   } else {
55648a46eb9SPierre Jolivet     if (sred->pctelescope_matcreate_type) PetscCall(sred->pctelescope_matcreate_type(pc, sred, MAT_REUSE_MATRIX, &sred->Bred));
5571e07b27eSBarry Smith   }
5581e07b27eSBarry Smith 
5591e07b27eSBarry Smith   /* common - no construction */
56057f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
5619566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(sred->ksp, sred->Bred, sred->Bred));
56248a46eb9SPierre Jolivet     if (pc->setfromoptionscalled && !pc->setupcalled) PetscCall(KSPSetFromOptions(sred->ksp));
5631e07b27eSBarry Smith   }
5643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5651e07b27eSBarry Smith }
5661e07b27eSBarry Smith 
PCApply_Telescope(PC pc,Vec x,Vec y)567d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Telescope(PC pc, Vec x, Vec y)
568d71ae5a4SJacob Faibussowitsch {
5691e07b27eSBarry Smith   PC_Telescope       sred = (PC_Telescope)pc->data;
5701e07b27eSBarry Smith   Vec                xtmp, xred, yred;
57113c30530SDave May   PetscInt           i, st, ed;
5721e07b27eSBarry Smith   VecScatter         scatter;
5731e07b27eSBarry Smith   PetscScalar       *array;
5741e07b27eSBarry Smith   const PetscScalar *x_array;
5751e07b27eSBarry Smith 
5761e07b27eSBarry Smith   PetscFunctionBegin;
5779566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(citation, &cited));
578bf00f589SPatrick Sanan 
5791e07b27eSBarry Smith   xtmp    = sred->xtmp;
5801e07b27eSBarry Smith   scatter = sred->scatter;
5811e07b27eSBarry Smith   xred    = sred->xred;
5821e07b27eSBarry Smith   yred    = sred->yred;
5831e07b27eSBarry Smith 
5841e07b27eSBarry Smith   /* pull in vector x->xtmp */
5859566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, x, xtmp, INSERT_VALUES, SCATTER_FORWARD));
5869566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, x, xtmp, INSERT_VALUES, SCATTER_FORWARD));
5871e07b27eSBarry Smith 
588bf00f589SPatrick Sanan   /* copy vector entries into xred */
5899566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(xtmp, &x_array));
5901e07b27eSBarry Smith   if (xred) {
5911e07b27eSBarry Smith     PetscScalar *LA_xred;
5929566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(xred, &st, &ed));
5939566063dSJacob Faibussowitsch     PetscCall(VecGetArray(xred, &LA_xred));
594ad540459SPierre Jolivet     for (i = 0; i < ed - st; i++) LA_xred[i] = x_array[i];
5959566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(xred, &LA_xred));
5961e07b27eSBarry Smith   }
5979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(xtmp, &x_array));
5981e07b27eSBarry Smith   /* solve */
59957f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
6009566063dSJacob Faibussowitsch     PetscCall(KSPSolve(sred->ksp, xred, yred));
6019566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(sred->ksp, pc, yred));
6021e07b27eSBarry Smith   }
6031e07b27eSBarry Smith   /* return vector */
6049566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xtmp, &array));
6051e07b27eSBarry Smith   if (yred) {
6061e07b27eSBarry Smith     const PetscScalar *LA_yred;
6079566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(yred, &st, &ed));
6089566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(yred, &LA_yred));
609ad540459SPierre Jolivet     for (i = 0; i < ed - st; i++) array[i] = LA_yred[i];
6109566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(yred, &LA_yred));
6111e07b27eSBarry Smith   }
6129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xtmp, &array));
6139566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scatter, xtmp, y, INSERT_VALUES, SCATTER_REVERSE));
6149566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scatter, xtmp, y, INSERT_VALUES, SCATTER_REVERSE));
6153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6161e07b27eSBarry Smith }
6171e07b27eSBarry Smith 
PCApplyRichardson_Telescope(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt its,PetscBool zeroguess,PetscInt * outits,PCRichardsonConvergedReason * reason)618d71ae5a4SJacob Faibussowitsch 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)
619d71ae5a4SJacob Faibussowitsch {
620f650675bSDave May   PC_Telescope       sred = (PC_Telescope)pc->data;
621a1d91a28SDave May   Vec                xtmp, yred;
622f650675bSDave May   PetscInt           i, st, ed;
623f650675bSDave May   VecScatter         scatter;
624f650675bSDave May   const PetscScalar *x_array;
625f650675bSDave May   PetscBool          default_init_guess_value;
626f650675bSDave May 
627f650675bSDave May   PetscFunctionBegin;
628f650675bSDave May   xtmp    = sred->xtmp;
629f650675bSDave May   scatter = sred->scatter;
630f650675bSDave May   yred    = sred->yred;
631f650675bSDave May 
63208401ef6SPierre Jolivet   PetscCheck(its <= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "PCApplyRichardson_Telescope only supports max_it = 1");
633f650675bSDave May   *reason = (PCRichardsonConvergedReason)0;
634f650675bSDave May 
635f650675bSDave May   if (!zeroguess) {
6369566063dSJacob Faibussowitsch     PetscCall(PetscInfo(pc, "PCTelescope: Scattering y for non-zero initial guess\n"));
637f650675bSDave May     /* pull in vector y->xtmp */
6389566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(scatter, y, xtmp, INSERT_VALUES, SCATTER_FORWARD));
6399566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(scatter, y, xtmp, INSERT_VALUES, SCATTER_FORWARD));
640f650675bSDave May 
641bf00f589SPatrick Sanan     /* copy vector entries into xred */
6429566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(xtmp, &x_array));
643f650675bSDave May     if (yred) {
644f650675bSDave May       PetscScalar *LA_yred;
6459566063dSJacob Faibussowitsch       PetscCall(VecGetOwnershipRange(yred, &st, &ed));
6469566063dSJacob Faibussowitsch       PetscCall(VecGetArray(yred, &LA_yred));
647ad540459SPierre Jolivet       for (i = 0; i < ed - st; i++) LA_yred[i] = x_array[i];
6489566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(yred, &LA_yred));
649f650675bSDave May     }
6509566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(xtmp, &x_array));
651f650675bSDave May   }
652f650675bSDave May 
65357f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
6549566063dSJacob Faibussowitsch     PetscCall(KSPGetInitialGuessNonzero(sred->ksp, &default_init_guess_value));
6559566063dSJacob Faibussowitsch     if (!zeroguess) PetscCall(KSPSetInitialGuessNonzero(sred->ksp, PETSC_TRUE));
656f650675bSDave May   }
657f650675bSDave May 
6589566063dSJacob Faibussowitsch   PetscCall(PCApply_Telescope(pc, x, y));
659f650675bSDave May 
66048a46eb9SPierre Jolivet   if (PCTelescope_isActiveRank(sred)) PetscCall(KSPSetInitialGuessNonzero(sred->ksp, default_init_guess_value));
661f650675bSDave May 
662f650675bSDave May   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
663f650675bSDave May   *outits = 1;
6643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
665f650675bSDave May }
666f650675bSDave May 
PCReset_Telescope(PC pc)667d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Telescope(PC pc)
668d71ae5a4SJacob Faibussowitsch {
6691e07b27eSBarry Smith   PC_Telescope sred = (PC_Telescope)pc->data;
6701e07b27eSBarry Smith 
671362febeeSStefano Zampini   PetscFunctionBegin;
6729566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&sred->isin));
6739566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&sred->scatter));
6749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sred->xred));
6759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sred->yred));
6769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&sred->xtmp));
6779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sred->Bred));
6789566063dSJacob Faibussowitsch   PetscCall(KSPReset(sred->ksp));
6791baa6e33SBarry Smith   if (sred->pctelescope_reset_type) PetscCall(sred->pctelescope_reset_type(pc));
6803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6811e07b27eSBarry Smith }
6821e07b27eSBarry Smith 
PCDestroy_Telescope(PC pc)683d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Telescope(PC pc)
684d71ae5a4SJacob Faibussowitsch {
6851e07b27eSBarry Smith   PC_Telescope sred = (PC_Telescope)pc->data;
6861e07b27eSBarry Smith 
6871e07b27eSBarry Smith   PetscFunctionBegin;
6889566063dSJacob Faibussowitsch   PetscCall(PCReset_Telescope(pc));
6899566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&sred->ksp));
6909566063dSJacob Faibussowitsch   PetscCall(PetscSubcommDestroy(&sred->psubcomm));
6919566063dSJacob Faibussowitsch   PetscCall(PetscFree(sred->dm_ctx));
6922e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetKSP_C", NULL));
6932e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetSubcommType_C", NULL));
6942e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetSubcommType_C", NULL));
6952e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetReductionFactor_C", NULL));
6962e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetReductionFactor_C", NULL));
6972e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetIgnoreDM_C", NULL));
6982e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetIgnoreDM_C", NULL));
6992e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetIgnoreKSPComputeOperators_C", NULL));
7002e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetIgnoreKSPComputeOperators_C", NULL));
7012e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetDM_C", NULL));
7022e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetUseCoarseDM_C", NULL));
7032e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetUseCoarseDM_C", NULL));
7049566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
7053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7061e07b27eSBarry Smith }
7071e07b27eSBarry Smith 
PCSetFromOptions_Telescope(PC pc,PetscOptionItems PetscOptionsObject)708ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_Telescope(PC pc, PetscOptionItems PetscOptionsObject)
709d71ae5a4SJacob Faibussowitsch {
7101e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
7111e07b27eSBarry Smith   MPI_Comm         comm;
7121e07b27eSBarry Smith   PetscMPIInt      size;
71348a10b22SPatrick Sanan   PetscBool        flg;
71448a10b22SPatrick Sanan   PetscSubcommType subcommtype;
7151e07b27eSBarry Smith 
7161e07b27eSBarry Smith   PetscFunctionBegin;
7179566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
7189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
719d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Telescope options");
7209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_telescope_subcomm_type", "Subcomm type (interlaced or contiguous)", "PCTelescopeSetSubcommType", PetscSubcommTypes, (PetscEnum)sred->subcommtype, (PetscEnum *)&subcommtype, &flg));
7211baa6e33SBarry Smith   if (flg) PetscCall(PCTelescopeSetSubcommType(pc, subcommtype));
7229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_telescope_reduction_factor", "Factor to reduce comm size by", "PCTelescopeSetReductionFactor", sred->redfactor, &sred->redfactor, NULL));
72308401ef6SPierre Jolivet   PetscCheck(sred->redfactor <= size, comm, PETSC_ERR_ARG_WRONG, "-pc_telescope_reduction_factor <= comm size");
7249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_telescope_ignore_dm", "Ignore any DM attached to the PC", "PCTelescopeSetIgnoreDM", sred->ignore_dm, &sred->ignore_dm, NULL));
7259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_telescope_ignore_kspcomputeoperators", "Ignore method used to compute A", "PCTelescopeSetIgnoreKSPComputeOperators", sred->ignore_kspcomputeoperators, &sred->ignore_kspcomputeoperators, NULL));
7269566063dSJacob 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));
727d0609cedSBarry Smith   PetscOptionsHeadEnd();
7283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7291e07b27eSBarry Smith }
7301e07b27eSBarry Smith 
731*5c7eeb11SPierre Jolivet /* PC implementation specific API's */
7321e07b27eSBarry Smith 
PCTelescopeGetKSP_Telescope(PC pc,KSP * ksp)733d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc, KSP *ksp)
734d71ae5a4SJacob Faibussowitsch {
7351e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
7364d86920dSPierre Jolivet 
737bd49479cSSatish Balay   PetscFunctionBegin;
7381e07b27eSBarry Smith   if (ksp) *ksp = red->ksp;
7393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7401e07b27eSBarry Smith }
7411e07b27eSBarry Smith 
PCTelescopeGetSubcommType_Telescope(PC pc,PetscSubcommType * subcommtype)742d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetSubcommType_Telescope(PC pc, PetscSubcommType *subcommtype)
743d71ae5a4SJacob Faibussowitsch {
74448a10b22SPatrick Sanan   PC_Telescope red = (PC_Telescope)pc->data;
7454d86920dSPierre Jolivet 
74648a10b22SPatrick Sanan   PetscFunctionBegin;
74748a10b22SPatrick Sanan   if (subcommtype) *subcommtype = red->subcommtype;
7483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74948a10b22SPatrick Sanan }
75048a10b22SPatrick Sanan 
PCTelescopeSetSubcommType_Telescope(PC pc,PetscSubcommType subcommtype)751d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeSetSubcommType_Telescope(PC pc, PetscSubcommType subcommtype)
752d71ae5a4SJacob Faibussowitsch {
75348a10b22SPatrick Sanan   PC_Telescope red = (PC_Telescope)pc->data;
75448a10b22SPatrick Sanan 
75548a10b22SPatrick Sanan   PetscFunctionBegin;
75628b400f6SJacob 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.");
75748a10b22SPatrick Sanan   red->subcommtype = subcommtype;
7583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75948a10b22SPatrick Sanan }
76048a10b22SPatrick Sanan 
PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt * fact)761d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc, PetscInt *fact)
762d71ae5a4SJacob Faibussowitsch {
7631e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
7644d86920dSPierre Jolivet 
765bd49479cSSatish Balay   PetscFunctionBegin;
7661e07b27eSBarry Smith   if (fact) *fact = red->redfactor;
7673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7681e07b27eSBarry Smith }
7691e07b27eSBarry Smith 
PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact)770d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc, PetscInt fact)
771d71ae5a4SJacob Faibussowitsch {
7721e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
7731e07b27eSBarry Smith   PetscMPIInt  size;
7741e07b27eSBarry Smith 
775bd49479cSSatish Balay   PetscFunctionBegin;
7769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
77763a3b9bcSJacob Faibussowitsch   PetscCheck(fact > 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Reduction factor of telescoping PC %" PetscInt_FMT " must be positive", fact);
77863a3b9bcSJacob Faibussowitsch   PetscCheck(fact <= size, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Reduction factor of telescoping PC %" PetscInt_FMT " must be <= comm.size", fact);
7791e07b27eSBarry Smith   red->redfactor = fact;
7803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7811e07b27eSBarry Smith }
7821e07b27eSBarry Smith 
PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool * v)783d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc, PetscBool *v)
784d71ae5a4SJacob Faibussowitsch {
7851e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
7864d86920dSPierre Jolivet 
787bd49479cSSatish Balay   PetscFunctionBegin;
7881e07b27eSBarry Smith   if (v) *v = red->ignore_dm;
7893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7901e07b27eSBarry Smith }
79148a10b22SPatrick Sanan 
PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v)792d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc, PetscBool v)
793d71ae5a4SJacob Faibussowitsch {
7941e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
7954d86920dSPierre Jolivet 
796bd49479cSSatish Balay   PetscFunctionBegin;
7971e07b27eSBarry Smith   red->ignore_dm = v;
7983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7991e07b27eSBarry Smith }
8001e07b27eSBarry Smith 
PCTelescopeGetUseCoarseDM_Telescope(PC pc,PetscBool * v)801d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetUseCoarseDM_Telescope(PC pc, PetscBool *v)
802d71ae5a4SJacob Faibussowitsch {
8038d9f7141SDave May   PC_Telescope red = (PC_Telescope)pc->data;
8044d86920dSPierre Jolivet 
8058d9f7141SDave May   PetscFunctionBegin;
8068d9f7141SDave May   if (v) *v = red->use_coarse_dm;
8073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8088d9f7141SDave May }
8098d9f7141SDave May 
PCTelescopeSetUseCoarseDM_Telescope(PC pc,PetscBool v)810d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeSetUseCoarseDM_Telescope(PC pc, PetscBool v)
811d71ae5a4SJacob Faibussowitsch {
8128d9f7141SDave May   PC_Telescope red = (PC_Telescope)pc->data;
8134d86920dSPierre Jolivet 
8148d9f7141SDave May   PetscFunctionBegin;
8158d9f7141SDave May   red->use_coarse_dm = v;
8163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8178d9f7141SDave May }
8188d9f7141SDave May 
PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool * v)819d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc, PetscBool *v)
820d71ae5a4SJacob Faibussowitsch {
8210ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
8224d86920dSPierre Jolivet 
8230ae7c45bSDave May   PetscFunctionBegin;
8240ae7c45bSDave May   if (v) *v = red->ignore_kspcomputeoperators;
8253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8260ae7c45bSDave May }
82748a10b22SPatrick Sanan 
PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool v)828d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc, PetscBool v)
829d71ae5a4SJacob Faibussowitsch {
8300ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
8314d86920dSPierre Jolivet 
8320ae7c45bSDave May   PetscFunctionBegin;
8330ae7c45bSDave May   red->ignore_kspcomputeoperators = v;
8343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8350ae7c45bSDave May }
8360ae7c45bSDave May 
PCTelescopeGetDM_Telescope(PC pc,DM * dm)837d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc, DM *dm)
838d71ae5a4SJacob Faibussowitsch {
8391e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
8404d86920dSPierre Jolivet 
841bd49479cSSatish Balay   PetscFunctionBegin;
8421e07b27eSBarry Smith   *dm = private_PCTelescopeGetSubDM(red);
8433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8441e07b27eSBarry Smith }
8451e07b27eSBarry Smith 
8461e07b27eSBarry Smith /*@
847f1580f4eSBarry Smith   PCTelescopeGetKSP - Gets the `KSP` created by the telescoping `PC`.
8481e07b27eSBarry Smith 
8491e07b27eSBarry Smith   Not Collective
8501e07b27eSBarry Smith 
8511e07b27eSBarry Smith   Input Parameter:
8521e07b27eSBarry Smith . pc - the preconditioner context
8531e07b27eSBarry Smith 
8541e07b27eSBarry Smith   Output Parameter:
85504c3f3b8SBarry Smith . subksp - the `KSP` defined on the smaller set of processes
8561e07b27eSBarry Smith 
8571e07b27eSBarry Smith   Level: advanced
8581e07b27eSBarry Smith 
859562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `KSP`, `PCTELESCOPE`
8601e07b27eSBarry Smith @*/
PCTelescopeGetKSP(PC pc,KSP * subksp)861d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetKSP(PC pc, KSP *subksp)
862d71ae5a4SJacob Faibussowitsch {
863bd49479cSSatish Balay   PetscFunctionBegin;
864cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetKSP_C", (PC, KSP *), (pc, subksp));
8653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8661e07b27eSBarry Smith }
8671e07b27eSBarry Smith 
8681e07b27eSBarry Smith /*@
86904c3f3b8SBarry Smith   PCTelescopeGetReductionFactor - Gets the factor by which the original number of MPI processes has been reduced by that was set by
87004c3f3b8SBarry Smith   `PCTelescopeSetReductionFactor()`
8711e07b27eSBarry Smith 
8721e07b27eSBarry Smith   Not Collective
8731e07b27eSBarry Smith 
8741e07b27eSBarry Smith   Input Parameter:
8751e07b27eSBarry Smith . pc - the preconditioner context
8761e07b27eSBarry Smith 
8771e07b27eSBarry Smith   Output Parameter:
8781e07b27eSBarry Smith . fact - the reduction factor
8791e07b27eSBarry Smith 
8801e07b27eSBarry Smith   Level: advanced
8811e07b27eSBarry Smith 
882562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCTELESCOPE`, `PCTelescopeSetReductionFactor()`
8831e07b27eSBarry Smith @*/
PCTelescopeGetReductionFactor(PC pc,PetscInt * fact)884d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetReductionFactor(PC pc, PetscInt *fact)
885d71ae5a4SJacob Faibussowitsch {
886bd49479cSSatish Balay   PetscFunctionBegin;
887cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetReductionFactor_C", (PC, PetscInt *), (pc, fact));
8883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8891e07b27eSBarry Smith }
8901e07b27eSBarry Smith 
8911e07b27eSBarry Smith /*@
8920b4b7b1cSBarry Smith   PCTelescopeSetReductionFactor - Sets the factor by which the original number of MPI processes will been reduced by when
8930b4b7b1cSBarry Smith   constructing the subcommunicator to be used with the `PCTELESCOPE`.
8941e07b27eSBarry Smith 
8951e07b27eSBarry Smith   Not Collective
8961e07b27eSBarry Smith 
8971e07b27eSBarry Smith   Input Parameter:
8981e07b27eSBarry Smith . pc - the preconditioner context
8991e07b27eSBarry Smith 
9001e07b27eSBarry Smith   Output Parameter:
9011e07b27eSBarry Smith . fact - the reduction factor
9021e07b27eSBarry Smith 
9031e07b27eSBarry Smith   Level: advanced
9041e07b27eSBarry Smith 
905562efe2eSBarry Smith .seealso: [](ch_ksp), `PCTELESCOPE`, `PCTelescopeGetReductionFactor()`
9061e07b27eSBarry Smith @*/
PCTelescopeSetReductionFactor(PC pc,PetscInt fact)907d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeSetReductionFactor(PC pc, PetscInt fact)
908d71ae5a4SJacob Faibussowitsch {
909bd49479cSSatish Balay   PetscFunctionBegin;
910cac4c232SBarry Smith   PetscTryMethod(pc, "PCTelescopeSetReductionFactor_C", (PC, PetscInt), (pc, fact));
9113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9121e07b27eSBarry Smith }
9131e07b27eSBarry Smith 
9141e07b27eSBarry Smith /*@
91504c3f3b8SBarry Smith   PCTelescopeGetIgnoreDM - Get the flag indicating if any `DM` attached to the `PC` will be used in constructing the `PC` on the
91604c3f3b8SBarry Smith   reduced number of MPI processes
9171e07b27eSBarry Smith 
9181e07b27eSBarry Smith   Not Collective
9191e07b27eSBarry Smith 
9201e07b27eSBarry Smith   Input Parameter:
9211e07b27eSBarry Smith . pc - the preconditioner context
9221e07b27eSBarry Smith 
9231e07b27eSBarry Smith   Output Parameter:
9241e07b27eSBarry Smith . v - the flag
9251e07b27eSBarry Smith 
9261e07b27eSBarry Smith   Level: advanced
9271e07b27eSBarry Smith 
928562efe2eSBarry Smith .seealso: [](ch_ksp), `DM`, `PCTELESCOPE`, `PCTelescopeSetIgnoreDM()`
9291e07b27eSBarry Smith @*/
PCTelescopeGetIgnoreDM(PC pc,PetscBool * v)930d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetIgnoreDM(PC pc, PetscBool *v)
931d71ae5a4SJacob Faibussowitsch {
932bd49479cSSatish Balay   PetscFunctionBegin;
933cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetIgnoreDM_C", (PC, PetscBool *), (pc, v));
9343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9351e07b27eSBarry Smith }
9361e07b27eSBarry Smith 
9371e07b27eSBarry Smith /*@
93804c3f3b8SBarry Smith   PCTelescopeSetIgnoreDM - Set a flag to ignore any `DM` attached to the `PC` when constructing the `PC` on the
93904c3f3b8SBarry Smith   reduced number of MPI processes
9401e07b27eSBarry Smith 
9411e07b27eSBarry Smith   Not Collective
9421e07b27eSBarry Smith 
9431e07b27eSBarry Smith   Input Parameter:
9441e07b27eSBarry Smith . pc - the preconditioner context
9451e07b27eSBarry Smith 
9461e07b27eSBarry Smith   Output Parameter:
94704c3f3b8SBarry Smith . v - Use `PETSC_TRUE` to ignore any `DM`
9481e07b27eSBarry Smith 
9491e07b27eSBarry Smith   Level: advanced
9501e07b27eSBarry Smith 
951562efe2eSBarry Smith .seealso: [](ch_ksp), `DM`, `PCTELESCOPE`, `PCTelescopeGetIgnoreDM()`
9521e07b27eSBarry Smith @*/
PCTelescopeSetIgnoreDM(PC pc,PetscBool v)953d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeSetIgnoreDM(PC pc, PetscBool v)
954d71ae5a4SJacob Faibussowitsch {
955bd49479cSSatish Balay   PetscFunctionBegin;
956cac4c232SBarry Smith   PetscTryMethod(pc, "PCTelescopeSetIgnoreDM_C", (PC, PetscBool), (pc, v));
9573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9581e07b27eSBarry Smith }
9591e07b27eSBarry Smith 
9601e07b27eSBarry Smith /*@
96104c3f3b8SBarry Smith   PCTelescopeGetUseCoarseDM - Get the flag indicating if the coarse `DM` attached to `DM` associated with the `PC` will be used in constructing
96204c3f3b8SBarry Smith   the `PC` on the reduced number of MPI processes
9638d9f7141SDave May 
9648d9f7141SDave May   Not Collective
9658d9f7141SDave May 
9668d9f7141SDave May   Input Parameter:
9678d9f7141SDave May . pc - the preconditioner context
9688d9f7141SDave May 
9698d9f7141SDave May   Output Parameter:
9708d9f7141SDave May . v - the flag
9718d9f7141SDave May 
9728d9f7141SDave May   Level: advanced
9738d9f7141SDave May 
974562efe2eSBarry Smith .seealso: [](ch_ksp), `DM`, `PCTELESCOPE`, `PCTelescopeSetIgnoreDM()`, `PCTelescopeSetUseCoarseDM()`
9758d9f7141SDave May @*/
PCTelescopeGetUseCoarseDM(PC pc,PetscBool * v)976d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetUseCoarseDM(PC pc, PetscBool *v)
977d71ae5a4SJacob Faibussowitsch {
9788d9f7141SDave May   PetscFunctionBegin;
979cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetUseCoarseDM_C", (PC, PetscBool *), (pc, v));
9803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9818d9f7141SDave May }
9828d9f7141SDave May 
9838d9f7141SDave May /*@
98404c3f3b8SBarry Smith   PCTelescopeSetUseCoarseDM - Set a flag to query the `DM` attached to the `PC` if it also has a coarse `DM` and utilize that `DM`
98504c3f3b8SBarry Smith   in constructing the `PC` on the reduced number of MPI processes
9868d9f7141SDave May 
9878d9f7141SDave May   Not Collective
9888d9f7141SDave May 
9898d9f7141SDave May   Input Parameter:
9908d9f7141SDave May . pc - the preconditioner context
9918d9f7141SDave May 
9928d9f7141SDave May   Output Parameter:
993f1580f4eSBarry Smith . v - Use `PETSC_FALSE` to ignore any coarse `DM`
9948d9f7141SDave May 
99504c3f3b8SBarry Smith   Level: advanced
99604c3f3b8SBarry Smith 
9978d9f7141SDave May   Notes:
99804c3f3b8SBarry Smith   When you have specified to use a coarse `DM`, the communicator used to create the sub-`KSP` within `PCTELESCOPE`
99904c3f3b8SBarry Smith   will be that of the coarse `DM`. Hence the flags `-pc_telescope_reduction_factor` and
10000b4b7b1cSBarry Smith   `-pc_telescope_subcomm_type` will not be used.
100104c3f3b8SBarry Smith 
1002f1580f4eSBarry Smith   It is required that the communicator associated with the parent (fine) and the coarse `DM` are of different sizes.
10030b4b7b1cSBarry Smith   An error will occur of the size if the communicator associated with the coarse `DM` is the same as that of the parent `DM`.
100404c3f3b8SBarry Smith   Furthermore, it is required that the communicator on the coarse `DM` is a sub-communicator of the parent.
10058d9f7141SDave May   This will be checked at the time the preconditioner is setup and an error will occur if
100604c3f3b8SBarry Smith   the coarse `DM` does not define a sub-communicator of that used by the parent `DM`.
10078d9f7141SDave May 
100804c3f3b8SBarry Smith   The particular Telescope setup invoked when using a coarse `DM` is agnostic with respect to the type of
1009f1580f4eSBarry Smith   the `DM` used (e.g. it supports `DMSHELL`, `DMPLEX`, etc).
10108d9f7141SDave May 
1011f1580f4eSBarry Smith   Support is currently only provided for the case when you are using `KSPSetComputeOperators()`
10128d9f7141SDave May 
101304c3f3b8SBarry Smith   The user is required to compose a function with the parent `DM` to facilitate the transfer of fields (`Vec`)
101404c3f3b8SBarry Smith   between the different decompositions defined by the fine and coarse `DM`s.
10158d9f7141SDave May   In the user code, this is achieved via
10168d9f7141SDave May .vb
10178d9f7141SDave May    {
10188d9f7141SDave May      DM dm_fine;
10198d9f7141SDave May      PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeFieldScatter",your_field_scatter_method);
10208d9f7141SDave May    }
10218d9f7141SDave May .ve
10228d9f7141SDave May   The signature of the user provided field scatter method is
10238d9f7141SDave May .vb
10248d9f7141SDave May    PetscErrorCode your_field_scatter_method(DM dm_fine,Vec x_fine,ScatterMode mode,DM dm_coarse,Vec x_coarse);
10258d9f7141SDave May .ve
102604c3f3b8SBarry Smith   The user must provide support for both mode `SCATTER_FORWARD` and mode `SCATTER_REVERSE`.
1027f1580f4eSBarry Smith   `SCATTER_FORWARD` implies the direction of transfer is from the parent (fine) `DM` to the coarse `DM`.
10288d9f7141SDave May 
102904c3f3b8SBarry Smith   Optionally, the user may also compose a function with the parent `DM` to facilitate the transfer
1030f1580f4eSBarry Smith   of state variables between the fine and coarse `DM`s.
10318d9f7141SDave May   In the context of a finite element discretization, an example state variable might be
10328d9f7141SDave May   values associated with quadrature points within each element.
10338d9f7141SDave May   A user provided state scatter method is composed via
10348d9f7141SDave May .vb
10358d9f7141SDave May    {
10368d9f7141SDave May      DM dm_fine;
10378d9f7141SDave May      PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeStateScatter",your_state_scatter_method);
10388d9f7141SDave May    }
10398d9f7141SDave May .ve
10408d9f7141SDave May   The signature of the user provided state scatter method is
10418d9f7141SDave May .vb
10428d9f7141SDave May    PetscErrorCode your_state_scatter_method(DM dm_fine,ScatterMode mode,DM dm_coarse);
10438d9f7141SDave May .ve
1044f1580f4eSBarry Smith   `SCATTER_FORWARD` implies the direction of transfer is from the fine `DM` to the coarse `DM`.
1045f1580f4eSBarry Smith   The user is only required to support mode = `SCATTER_FORWARD`.
10468d9f7141SDave May   No assumption is made about the data type of the state variables.
1047f1580f4eSBarry Smith   These must be managed by the user and must be accessible from the `DM`.
10488d9f7141SDave May 
1049f1580f4eSBarry Smith   Care must be taken in defining the user context passed to `KSPSetComputeOperators()` which is to be
1050f1580f4eSBarry Smith   associated with the sub-`KSP` residing within `PCTELESCOPE`.
1051f1580f4eSBarry Smith   In general, `PCTELESCOPE` assumes that the context on the fine and coarse `DM` used with
1052f1580f4eSBarry Smith   `KSPSetComputeOperators()` should be "similar" in type or origin.
1053f1580f4eSBarry Smith   Specifically the following rules are used to infer what context on the sub-`KSP` should be.
10548d9f7141SDave May 
1055f1580f4eSBarry Smith   First the contexts from the `KSP` and the fine and coarse `DM`s are retrieved.
1056f1580f4eSBarry Smith   Note that the special case of a `DMSHELL` context is queried.
10578d9f7141SDave May 
10588d9f7141SDave May .vb
10598d9f7141SDave May    DMKSPGetComputeOperators(dm_fine,&dmfine_kspfunc,&dmfine_kspctx);
10608d9f7141SDave May    DMGetApplicationContext(dm_fine,&dmfine_appctx);
10618d9f7141SDave May    DMShellGetContext(dm_fine,&dmfine_shellctx);
10628d9f7141SDave May 
10638d9f7141SDave May    DMGetApplicationContext(dm_coarse,&dmcoarse_appctx);
10648d9f7141SDave May    DMShellGetContext(dm_coarse,&dmcoarse_shellctx);
10658d9f7141SDave May .ve
10668d9f7141SDave May 
106704c3f3b8SBarry Smith   The following rules are then enforced\:
10688d9f7141SDave May 
106904c3f3b8SBarry Smith   1. If `dmfine_kspctx` = `NULL`, then we provide a `NULL` pointer as the context for the sub-`KSP`\:
107004c3f3b8SBarry Smith   `KSPSetComputeOperators`(`sub_ksp`,`dmfine_kspfunc`,`NULL`);
10718d9f7141SDave May 
107204c3f3b8SBarry Smith   2. If `dmfine_kspctx` != `NULL` and `dmfine_kspctx` == `dmfine_appctx`,
1073feefa0e1SJacob Faibussowitsch 
107404c3f3b8SBarry Smith   check that `dmcoarse_appctx` is also non-`NULL`. If this is true, then\:
107504c3f3b8SBarry Smith   `KSPSetComputeOperators`(`sub_ksp`,`dmfine_kspfunc`,`dmcoarse_appctx`);
10768d9f7141SDave May 
107704c3f3b8SBarry Smith   3. If `dmfine_kspctx` != `NULL` and `dmfine_kspctx` == `dmfine_shellctx`,
1078feefa0e1SJacob Faibussowitsch 
107904c3f3b8SBarry Smith   check that `dmcoarse_shellctx` is also non-`NULL`. If this is true, then\:
108004c3f3b8SBarry Smith   `KSPSetComputeOperators`(`sub_ksp`,`dmfine_kspfunc`,`dmcoarse_shellctx`);
10818d9f7141SDave May 
1082f1580f4eSBarry Smith   If neither of the above three tests passed, then `PCTELESCOPE` cannot safely determine what
1083f1580f4eSBarry Smith   context should be provided to `KSPSetComputeOperators()` for use with the sub-`KSP`.
10848d9f7141SDave May   In this case, an additional mechanism is provided via a composed function which will return
10858d9f7141SDave May   the actual context to be used. To use this feature you must compose the "getter" function
1086f1580f4eSBarry Smith   with the coarse `DM`, e.g.
10878d9f7141SDave May .vb
10888d9f7141SDave May    {
10898d9f7141SDave May      DM dm_coarse;
10908d9f7141SDave May      PetscObjectCompose((PetscObject)dm_coarse,"PCTelescopeGetCoarseDMKSPContext",your_coarse_context_getter);
10918d9f7141SDave May    }
10928d9f7141SDave May .ve
10938d9f7141SDave May   The signature of the user provided method is
10948d9f7141SDave May .vb
10958d9f7141SDave May    PetscErrorCode your_coarse_context_getter(DM dm_coarse,void **your_kspcontext);
10968d9f7141SDave May .ve
10978d9f7141SDave May 
1098562efe2eSBarry Smith .seealso: [](ch_ksp), `DM`, `PCTELESCOPE`, `PCTelescopeSetIgnoreDM()`
10998d9f7141SDave May @*/
PCTelescopeSetUseCoarseDM(PC pc,PetscBool v)1100d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeSetUseCoarseDM(PC pc, PetscBool v)
1101d71ae5a4SJacob Faibussowitsch {
11028d9f7141SDave May   PetscFunctionBegin;
1103cac4c232SBarry Smith   PetscTryMethod(pc, "PCTelescopeSetUseCoarseDM_C", (PC, PetscBool), (pc, v));
11043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11058d9f7141SDave May }
11068d9f7141SDave May 
11078d9f7141SDave May /*@
110804c3f3b8SBarry Smith   PCTelescopeGetIgnoreKSPComputeOperators - Get the flag indicating if `KSPComputeOperators()` will be used to construct
110904c3f3b8SBarry Smith   the matrix on the reduced number of MPI processes
11100ae7c45bSDave May 
11110ae7c45bSDave May   Not Collective
11120ae7c45bSDave May 
11130ae7c45bSDave May   Input Parameter:
11140ae7c45bSDave May . pc - the preconditioner context
11150ae7c45bSDave May 
11160ae7c45bSDave May   Output Parameter:
11170ae7c45bSDave May . v - the flag
11180ae7c45bSDave May 
11190ae7c45bSDave May   Level: advanced
11200ae7c45bSDave May 
1121562efe2eSBarry Smith .seealso: [](ch_ksp), `PCTELESCOPE`, `PCTelescopeSetIgnoreDM()`, `PCTelescopeSetUseCoarseDM()`, `PCTelescopeSetIgnoreKSPComputeOperators()`
11220ae7c45bSDave May @*/
PCTelescopeGetIgnoreKSPComputeOperators(PC pc,PetscBool * v)1123d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC pc, PetscBool *v)
1124d71ae5a4SJacob Faibussowitsch {
11250ae7c45bSDave May   PetscFunctionBegin;
1126cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetIgnoreKSPComputeOperators_C", (PC, PetscBool *), (pc, v));
11273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11280ae7c45bSDave May }
11290ae7c45bSDave May 
11300ae7c45bSDave May /*@
113104c3f3b8SBarry Smith   PCTelescopeSetIgnoreKSPComputeOperators - Set a flag to have `PCTELESCOPE` ignore the function provided to `KSPComputeOperators()` in
113204c3f3b8SBarry Smith   constructint the matrix on the reduced number of MPI processes
11330ae7c45bSDave May 
11340ae7c45bSDave May   Not Collective
11350ae7c45bSDave May 
11360ae7c45bSDave May   Input Parameter:
11370ae7c45bSDave May . pc - the preconditioner context
11380ae7c45bSDave May 
11390ae7c45bSDave May   Output Parameter:
114004c3f3b8SBarry Smith . v - Use `PETSC_TRUE` to ignore the function (if defined) set via `KSPSetComputeOperators()` on `pc`
11410ae7c45bSDave May 
11420ae7c45bSDave May   Level: advanced
11430ae7c45bSDave May 
1144562efe2eSBarry Smith .seealso: [](ch_ksp), `PCTELESCOPE`, `PCTelescopeSetIgnoreDM()`, `PCTelescopeSetUseCoarseDM()`, `PCTelescopeGetIgnoreKSPComputeOperators()`
11450ae7c45bSDave May @*/
PCTelescopeSetIgnoreKSPComputeOperators(PC pc,PetscBool v)1146d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC pc, PetscBool v)
1147d71ae5a4SJacob Faibussowitsch {
11480ae7c45bSDave May   PetscFunctionBegin;
1149cac4c232SBarry Smith   PetscTryMethod(pc, "PCTelescopeSetIgnoreKSPComputeOperators_C", (PC, PetscBool), (pc, v));
11503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11510ae7c45bSDave May }
11520ae7c45bSDave May 
11530ae7c45bSDave May /*@
1154f1580f4eSBarry Smith   PCTelescopeGetDM - Get the re-partitioned `DM` attached to the sub-`KSP`.
11551e07b27eSBarry Smith 
11561e07b27eSBarry Smith   Not Collective
11571e07b27eSBarry Smith 
11581e07b27eSBarry Smith   Input Parameter:
11591e07b27eSBarry Smith . pc - the preconditioner context
11601e07b27eSBarry Smith 
11611e07b27eSBarry Smith   Output Parameter:
116204c3f3b8SBarry Smith . subdm - The re-partitioned `DM`
11631e07b27eSBarry Smith 
11641e07b27eSBarry Smith   Level: advanced
11651e07b27eSBarry Smith 
1166562efe2eSBarry Smith .seealso: [](ch_ksp), `DM`, `PCTELESCOPE`, `PCTelescopeSetIgnoreDM()`, `PCTelescopeSetUseCoarseDM()`, `PCTelescopeGetIgnoreKSPComputeOperators()`
11671e07b27eSBarry Smith @*/
PCTelescopeGetDM(PC pc,DM * subdm)1168d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetDM(PC pc, DM *subdm)
1169d71ae5a4SJacob Faibussowitsch {
1170bd49479cSSatish Balay   PetscFunctionBegin;
1171cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetDM_C", (PC, DM *), (pc, subdm));
11723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11731e07b27eSBarry Smith }
11741e07b27eSBarry Smith 
117548a10b22SPatrick Sanan /*@
11760b4b7b1cSBarry Smith   PCTelescopeSetSubcommType - set subcommunicator type `PetscSubcommType` (interlaced or contiguous) to be used when
11770b4b7b1cSBarry Smith   the subcommunicator is generated from the given `PC`
117848a10b22SPatrick Sanan 
117948a10b22SPatrick Sanan   Logically Collective
118048a10b22SPatrick Sanan 
1181d8d19677SJose E. Roman   Input Parameters:
11821dae98e4SBarry Smith + pc          - the preconditioner context
1183f1580f4eSBarry Smith - subcommtype - the subcommunicator type (see `PetscSubcommType`)
118448a10b22SPatrick Sanan 
118548a10b22SPatrick Sanan   Level: advanced
118648a10b22SPatrick Sanan 
1187562efe2eSBarry Smith .seealso: [](ch_ksp), `PetscSubcommType`, `PetscSubcomm`, `PCTELESCOPE`, `PCTelescopeGetSubcommType()`
118848a10b22SPatrick Sanan @*/
PCTelescopeSetSubcommType(PC pc,PetscSubcommType subcommtype)1189d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeSetSubcommType(PC pc, PetscSubcommType subcommtype)
1190d71ae5a4SJacob Faibussowitsch {
119148a10b22SPatrick Sanan   PetscFunctionBegin;
1192cac4c232SBarry Smith   PetscTryMethod(pc, "PCTelescopeSetSubcommType_C", (PC, PetscSubcommType), (pc, subcommtype));
11933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119448a10b22SPatrick Sanan }
119548a10b22SPatrick Sanan 
119648a10b22SPatrick Sanan /*@
11970b4b7b1cSBarry Smith   PCTelescopeGetSubcommType - Get the subcommunicator type `PetscSubcommType` (interlaced or contiguous) set with `PCTelescopeSetSubcommType()`
119848a10b22SPatrick Sanan 
119948a10b22SPatrick Sanan   Not Collective
120048a10b22SPatrick Sanan 
120148a10b22SPatrick Sanan   Input Parameter:
120248a10b22SPatrick Sanan . pc - the preconditioner context
120348a10b22SPatrick Sanan 
120448a10b22SPatrick Sanan   Output Parameter:
1205f1580f4eSBarry Smith . subcommtype - the subcommunicator type (see `PetscSubcommType`)
120648a10b22SPatrick Sanan 
120748a10b22SPatrick Sanan   Level: advanced
120848a10b22SPatrick Sanan 
1209562efe2eSBarry Smith .seealso: [](ch_ksp), `PetscSubcomm`, `PetscSubcommType`, `PCTELESCOPE`, `PCTelescopeSetSubcommType()`
121048a10b22SPatrick Sanan @*/
PCTelescopeGetSubcommType(PC pc,PetscSubcommType * subcommtype)1211d71ae5a4SJacob Faibussowitsch PetscErrorCode PCTelescopeGetSubcommType(PC pc, PetscSubcommType *subcommtype)
1212d71ae5a4SJacob Faibussowitsch {
121348a10b22SPatrick Sanan   PetscFunctionBegin;
1214cac4c232SBarry Smith   PetscUseMethod(pc, "PCTelescopeGetSubcommType_C", (PC, PetscSubcommType *), (pc, subcommtype));
12153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
121648a10b22SPatrick Sanan }
121748a10b22SPatrick Sanan 
12181e07b27eSBarry Smith /*MC
12190b4b7b1cSBarry Smith    PCTELESCOPE - Runs a `KSP` solver on a sub-communicator {cite}`maysananruppknepleysmith2016` of the communicator used by the original `KSP`.
12200b4b7b1cSBarry Smith                  MPI processes not in the sub-communicator are idle during the solve. Usually used to solve the smaller coarser grid problems in multigrid
12210b4b7b1cSBarry Smith                  (`PCMG`) that could not be efficiently solved on the entire communication
12221e07b27eSBarry Smith 
1223f1580f4eSBarry Smith    Options Database Keys:
122400fea0ebSDave 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.
122504c3f3b8SBarry Smith .  -pc_telescope_ignore_dm                            - flag to indicate whether an attached `DM` should be ignored in constructing the new `PC`
122604c3f3b8SBarry Smith .  -pc_telescope_subcomm_type <interlaced,contiguous> - defines the selection of MPI processes on the sub-communicator. see `PetscSubcomm` for more information.
122704c3f3b8SBarry Smith .  -pc_telescope_ignore_kspcomputeoperators           - flag to indicate whether `KSPSetComputeOperators()` should be used on the sub-`KSP`.
1228f1580f4eSBarry Smith -  -pc_telescope_use_coarse_dm                        - flag to indicate whether the coarse `DM` should be used to define the sub-communicator.
12291e07b27eSBarry Smith 
12301e07b27eSBarry Smith    Level: advanced
12311e07b27eSBarry Smith 
12321e07b27eSBarry Smith    Notes:
1233f1580f4eSBarry Smith    Assuming that the parent preconditioner `PC` is defined on a communicator c, this implementation
123404c3f3b8SBarry Smith    creates a child sub-communicator (c') containing fewer MPI processes than the original parent preconditioner `PC`.
1235f1580f4eSBarry Smith    The preconditioner is deemed telescopic as it only calls `KSPSolve()` on a single
1236f1580f4eSBarry Smith    sub-communicator, in contrast with `PCREDUNDANT` which calls `KSPSolve()` on N sub-communicators.
123704c3f3b8SBarry Smith    This means there will be MPI processes which will be idle during the application of this preconditioner.
123804c3f3b8SBarry Smith    Additionally, in comparison with `PCREDUNDANT`, `PCTELESCOPE` can utilize an attached `DM` to construct `DM` dependent preconditioner, such as `PCMG`
12396fc41876SBarry Smith 
12400b4b7b1cSBarry Smith    The default type `KSPType` of the sub `KSP` (the `KSP` defined on c') is `KSPPREONLY`.
124100fea0ebSDave May 
1242f1580f4eSBarry Smith    There are three setup mechanisms for `PCTELESCOPE`. Features support by each type are described below.
12430b4b7b1cSBarry Smith    In the following, we will refer to the operators B and B', these are the `Bmat` provided to the `KSP` on the
124400fea0ebSDave May    communicators c and c' respectively.
124500fea0ebSDave May 
124600fea0ebSDave May    [1] Default setup
1247f1580f4eSBarry Smith    The sub-communicator c' is created via `PetscSubcommCreate()`.
12480b4b7b1cSBarry Smith    Any explicitly defined nullspace and near nullspace vectors attached to B with `MatSetNullSpace()` and `MatSetNearNullSpace()` are transferred to B'.
12490b4b7b1cSBarry Smith    Currently there is no support for nullspaces provided with `MatNullSpaceSetFunction()`).
1250f1580f4eSBarry Smith    No support is provided for `KSPSetComputeOperators()`.
125104c3f3b8SBarry Smith    Currently there is no support for the flag `-pc_use_amat`.
125200fea0ebSDave May 
1253f1580f4eSBarry Smith    [2] `DM` aware setup
12540b4b7b1cSBarry Smith    The sub-communicator c' is created via `PetscSubcommCreate()`.
1255f1580f4eSBarry Smith    If a `DM` is attached to the `PC`, it is re-partitioned on the sub-communicator c'.
12560b4b7b1cSBarry Smith    Both the `Bmat` operator and the right-hand side vector are permuted into the new DOF ordering defined by the re-partitioned `DM`.
1257f1580f4eSBarry Smith    Currently only support for re-partitioning a `DMDA` is provided.
12580b4b7b1cSBarry Smith    Any explicitly defined nullspace or near nullspace vectors attached to the original B with `MatSetNullSpace()`
12590b4b7b1cSBarry Smith    and `MatSetNearNullSpace()` are extracted, re-partitioned and set on B'
12600b4b7b1cSBarry Smith    (currently there is no support for nullspaces provided with `MatNullSpaceSetFunction()`).
1261f1580f4eSBarry Smith    Support is provided for `KSPSetComputeOperators()`. The user provided function and context is propagated to the sub `KSP`.
126200fea0ebSDave May    This is fragile since the user must ensure that their user context is valid for use on c'.
126304c3f3b8SBarry Smith    Currently there is no support for the flag `-pc_use_amat`.
12641e07b27eSBarry Smith 
1265f1580f4eSBarry Smith    [3] Coarse `DM` setup
1266f1580f4eSBarry Smith    If a `DM` (dmfine) is attached to the `PC`, dmfine is queried for a "coarse" `DM` (call this dmcoarse) via `DMGetCoarseDM()`.
1267f1580f4eSBarry Smith    `PCTELESCOPE` will interpret the coarse `DM` as being defined on a sub-communicator of c.
1268f1580f4eSBarry Smith    The communicator associated with dmcoarse will define the c' to be used within `PCTELESCOPE`.
1269f1580f4eSBarry Smith    `PCTELESCOPE` will check that c' is in fact a sub-communicator of c. If it is not, an error will be reported.
1270f1580f4eSBarry Smith    The intention of this setup type is that `PCTELESCOPE` will use an existing (e.g. user defined) communicator hierarchy, say as would be
127100fea0ebSDave May    available with using multi-grid on unstructured meshes.
127204c3f3b8SBarry Smith    This setup will not use the command line options `-pc_telescope_reduction_factor` or `-pc_telescope_subcomm_type`.
12730b4b7b1cSBarry Smith    Any explicitly defined nullspace or near nullspace vectors attached to the B are extracted, scattered into the correct ordering consistent
12740b4b7b1cSBarry Smith    with dmcoarse and set on B'
12750b4b7b1cSBarry Smith    (currently there is no support for nullspaces provided with `MatNullSpaceSetFunction()`).
1276f1580f4eSBarry Smith    There is no general method to permute field orderings, hence only `KSPSetComputeOperators()` is supported.
1277f1580f4eSBarry Smith    The user must use `PetscObjectComposeFunction()` with dmfine to define the method to scatter fields from dmfine to dmcoarse.
12780b4b7b1cSBarry Smith    Propagation of the user context for `KSPSetComputeOperators()` on the sub `KSP` is attempted by querying the `DM` contexts associated with
12790b4b7b1cSBarry Smith    dmfine and dmcoarse. Alternatively, the user may use `PetscObjectComposeFunction()` with dmcoarse to define a method which will return the appropriate user context for `KSPSetComputeOperators()`.
128004c3f3b8SBarry Smith    Currently there is no support for the flag `-pc_use_amat`.
128104c3f3b8SBarry Smith    This setup can be invoked by the option `-pc_telescope_use_coarse_dm` or by calling `PCTelescopeSetUseCoarseDM`(pc,`PETSC_TRUE`);
1282f1580f4eSBarry Smith    Further information about the user-provided methods required by this setup type are described here `PCTelescopeSetUseCoarseDM()`.
12836fc41876SBarry Smith 
12846fc41876SBarry Smith    Developer Notes:
128554c05997SPierre Jolivet    During `PCSetUp()`, the B operator is scattered onto c'.
1286f1580f4eSBarry Smith    Within `PCApply()`, the RHS vector (x) is scattered into a redundant vector, xred (defined on c').
1287f1580f4eSBarry Smith    Then, `KSPSolve()` is executed on the c' communicator.
12886fc41876SBarry Smith 
1289f1580f4eSBarry Smith    The communicator used within the telescoping preconditioner is defined by a `PetscSubcomm` using the INTERLACED
129004c3f3b8SBarry Smith    creation routine by default (this can be changed with `-pc_telescope_subcomm_type`). We run the sub `KSP` on only
129104c3f3b8SBarry Smith    the ranks within the communicator which have a color equal to zero.
12926fc41876SBarry Smith 
1293005d9f20SPatrick Sanan    The telescoping preconditioner is aware of nullspaces and near nullspaces which are attached to the B operator.
12948439623fSPatrick Sanan    In the case where B has a (near) nullspace attached, the (near) nullspace vectors are extracted from B and mapped into
1295005d9f20SPatrick Sanan    a new (near) nullspace, defined on the sub-communicator, which is attached to B' (the B operator which was scattered to c')
12966fc41876SBarry Smith 
129704c3f3b8SBarry Smith    The telescoping preconditioner can re-partition an attached `DM` if it is a `DMDA` (2D or 3D -
1298f1580f4eSBarry Smith    support for 1D `DMDA`s is not provided). If a `DMDA` is found, a topologically equivalent `DMDA` is created on c'
1299f1580f4eSBarry Smith    and this new `DM` is attached the sub `KSP`. The design of telescope is such that it should be possible to extend support
130004c3f3b8SBarry Smith    for re-partitioning other to `DM`'s (e.g. `DMPLEX`). The user can supply a flag to ignore attached DMs.
130104c3f3b8SBarry Smith    Alternatively, user-provided re-partitioned `DM`s can be used via `-pc_telescope_use_coarse_dm`.
13026fc41876SBarry Smith 
130304c3f3b8SBarry Smith    With the default setup mode, B' is defined by fusing rows (in order) associated with MPI processes common to c and c'.
13046fc41876SBarry Smith 
1305f1580f4eSBarry Smith    When a `DMDA` is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permutation of B
1306f1580f4eSBarry Smith    into the ordering defined by the `DMDA` on c', (ii) extracting the local chunks via `MatCreateSubMatrices()`, (iii) fusing the
1307f1580f4eSBarry Smith    locally (sequential) matrices defined on the ranks common to c and c' into B' using `MatCreateMPIMatConcatenateSeqMat()`
13086fc41876SBarry Smith 
13098439623fSPatrick Sanan    Limitations/improvements include the following.
1310f1580f4eSBarry Smith    `VecPlaceArray()` could be used within `PCApply()` to improve efficiency and reduce memory usage.
1311f1580f4eSBarry Smith    A unified mechanism to query for user contexts as required by `KSPSetComputeOperators()` and `MatNullSpaceSetFunction()`.
13126fc41876SBarry Smith 
1313f1580f4eSBarry Smith    The symmetric permutation used when a `DMDA` is encountered is performed via explicitly assembling a permutation matrix P,
1314f1580f4eSBarry Smith    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
1315f1580f4eSBarry Smith    `VecPermute()` does not support the use case required here. By computing P, one can permute both the operator and RHS in a
13166fc41876SBarry Smith    consistent manner.
13176fc41876SBarry Smith 
131800fea0ebSDave May    Mapping of vectors (default setup mode) is performed in the following way.
131900fea0ebSDave 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.
13208439623fSPatrick Sanan    Using the interlaced creation routine, the ranks in c with color = 0 will be rank 0 and 2.
132100fea0ebSDave May    We perform the scatter to the sub-communicator in the following way.
1322514db83aSDave May    [1] Given a vector x defined on communicator c
13236fc41876SBarry Smith 
1324514db83aSDave May .vb
1325514db83aSDave May    rank(c)  local values of x
1326514db83aSDave May    ------- ----------------------------------------
1327514db83aSDave May         0   [  0.0,  1.0,  2.0,  3.0,  4.0,  5.0 ]
1328514db83aSDave May         1   [  6.0,  7.0,  8.0,  9.0, 10.0, 11.0 ]
1329514db83aSDave May         2   [ 12.0, 13.0, 14.0, 15.0, 16.0, 17.0 ]
1330514db83aSDave May         3   [ 18.0, 19.0, 20.0, 21.0, 22.0, 23.0 ]
1331514db83aSDave May .ve
13326fc41876SBarry Smith 
1333514db83aSDave May    scatter into xtmp defined also on comm c, so that we have the following values
13346fc41876SBarry Smith 
1335514db83aSDave May .vb
1336514db83aSDave May    rank(c)  local values of xtmp
1337514db83aSDave May    ------- ----------------------------------------------------------------------------
1338514db83aSDave 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 ]
1339514db83aSDave May         1   [ ]
1340514db83aSDave 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 ]
1341514db83aSDave May         3   [ ]
1342514db83aSDave May .ve
13436fc41876SBarry Smith 
13446fc41876SBarry Smith    The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') have no values
13456fc41876SBarry Smith 
1346514db83aSDave 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'.
13476fc41876SBarry Smith    Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0.
13486fc41876SBarry Smith 
1349514db83aSDave May .vb
1350514db83aSDave May    rank(c')  local values of xred
1351514db83aSDave May    -------- ----------------------------------------------------------------------------
1352514db83aSDave 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 ]
1353514db83aSDave 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 ]
1354514db83aSDave May .ve
13551e07b27eSBarry Smith 
13561d27aa22SBarry Smith   Contributed by:
13571d27aa22SBarry Smith   Dave May
1358bf00f589SPatrick Sanan 
1359562efe2eSBarry Smith .seealso: [](ch_ksp), `PCTelescopeGetKSP()`, `PCTelescopeGetDM()`, `PCTelescopeGetReductionFactor()`, `PCTelescopeSetReductionFactor()`, `PCTelescopeGetIgnoreDM()`, `PCTelescopeSetIgnoreDM()`, `PCREDUNDANT`
13601e07b27eSBarry Smith M*/
PCCreate_Telescope(PC pc)1361d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc)
1362d71ae5a4SJacob Faibussowitsch {
13631e07b27eSBarry Smith   struct _PC_Telescope *sred;
13641e07b27eSBarry Smith 
13651e07b27eSBarry Smith   PetscFunctionBegin;
13664dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&sred));
13672a22aa42SDave May   sred->psubcomm                   = NULL;
136848a10b22SPatrick Sanan   sred->subcommtype                = PETSC_SUBCOMM_INTERLACED;
13692a22aa42SDave May   sred->subcomm                    = MPI_COMM_NULL;
13701e07b27eSBarry Smith   sred->redfactor                  = 1;
13711e07b27eSBarry Smith   sred->ignore_dm                  = PETSC_FALSE;
13727c5279cbSDave May   sred->ignore_kspcomputeoperators = PETSC_FALSE;
13738d9f7141SDave May   sred->use_coarse_dm              = PETSC_FALSE;
13741e07b27eSBarry Smith   pc->data                         = (void *)sred;
13751e07b27eSBarry Smith 
13761e07b27eSBarry Smith   pc->ops->apply           = PCApply_Telescope;
13771e07b27eSBarry Smith   pc->ops->applytranspose  = NULL;
1378f650675bSDave May   pc->ops->applyrichardson = PCApplyRichardson_Telescope;
13791e07b27eSBarry Smith   pc->ops->setup           = PCSetUp_Telescope;
13801e07b27eSBarry Smith   pc->ops->destroy         = PCDestroy_Telescope;
13811e07b27eSBarry Smith   pc->ops->reset           = PCReset_Telescope;
13821e07b27eSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Telescope;
13831e07b27eSBarry Smith   pc->ops->view            = PCView_Telescope;
13841e07b27eSBarry Smith 
13851e07b27eSBarry Smith   sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
13861e07b27eSBarry Smith   sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
13871e07b27eSBarry Smith   sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
13881e07b27eSBarry Smith   sred->pctelescope_reset_type              = NULL;
13891e07b27eSBarry Smith 
13909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetKSP_C", PCTelescopeGetKSP_Telescope));
13919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetSubcommType_C", PCTelescopeGetSubcommType_Telescope));
13929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetSubcommType_C", PCTelescopeSetSubcommType_Telescope));
13939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetReductionFactor_C", PCTelescopeGetReductionFactor_Telescope));
13949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetReductionFactor_C", PCTelescopeSetReductionFactor_Telescope));
13959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetIgnoreDM_C", PCTelescopeGetIgnoreDM_Telescope));
13969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetIgnoreDM_C", PCTelescopeSetIgnoreDM_Telescope));
13979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetIgnoreKSPComputeOperators_C", PCTelescopeGetIgnoreKSPComputeOperators_Telescope));
13989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetIgnoreKSPComputeOperators_C", PCTelescopeSetIgnoreKSPComputeOperators_Telescope));
13999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetDM_C", PCTelescopeGetDM_Telescope));
14009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeGetUseCoarseDM_C", PCTelescopeGetUseCoarseDM_Telescope));
14019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCTelescopeSetUseCoarseDM_C", PCTelescopeSetUseCoarseDM_Telescope));
14023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14031e07b27eSBarry Smith }
1404