xref: /petsc/src/ksp/pc/impls/telescope/telescope.c (revision a8d69d7b0124b1e6ce75950a93e6ff079980e86f)
11e07b27eSBarry Smith 
28d9f7141SDave May #include <petsc/private/petscimpl.h>
3120bdd93SDave May #include <petsc/private/matimpl.h>
46fc41876SBarry Smith #include <petsc/private/pcimpl.h>
51e07b27eSBarry Smith #include <petscksp.h> /*I "petscksp.h" I*/
61e07b27eSBarry Smith #include <petscdm.h> /*I "petscdm.h" I*/
7575a0592SBarry Smith #include "../src/ksp/pc/impls/telescope/telescope.h"
81e07b27eSBarry Smith 
9bf00f589SPatrick Sanan static PetscBool  cited = PETSC_FALSE;
10bf00f589SPatrick Sanan static const char citation[] =
11bf00f589SPatrick Sanan "@inproceedings{MaySananRuppKnepleySmith2016,\n"
12bf00f589SPatrick Sanan "  title     = {Extreme-Scale Multigrid Components within PETSc},\n"
13bf00f589SPatrick Sanan "  author    = {Dave A. May and Patrick Sanan and Karl Rupp and Matthew G. Knepley and Barry F. Smith},\n"
14bf00f589SPatrick Sanan "  booktitle = {Proceedings of the Platform for Advanced Scientific Computing Conference},\n"
15bf00f589SPatrick Sanan "  series    = {PASC '16},\n"
16bf00f589SPatrick Sanan "  isbn      = {978-1-4503-4126-4},\n"
17bf00f589SPatrick Sanan "  location  = {Lausanne, Switzerland},\n"
18bf00f589SPatrick Sanan "  pages     = {5:1--5:12},\n"
19bf00f589SPatrick Sanan "  articleno = {5},\n"
20bf00f589SPatrick Sanan "  numpages  = {12},\n"
21*a8d69d7bSBarry Smith "  url       = {https://doi.acm.org/10.1145/2929908.2929913},\n"
22bf00f589SPatrick Sanan "  doi       = {10.1145/2929908.2929913},\n"
23bf00f589SPatrick Sanan "  acmid     = {2929913},\n"
24bf00f589SPatrick Sanan "  publisher = {ACM},\n"
25bf00f589SPatrick Sanan "  address   = {New York, NY, USA},\n"
26bf00f589SPatrick Sanan "  keywords  = {GPU, HPC, agglomeration, coarse-level solver, multigrid, parallel computing, preconditioning},\n"
27bf00f589SPatrick Sanan "  year      = {2016}\n"
28bf00f589SPatrick Sanan "}\n";
29bf00f589SPatrick Sanan 
301e07b27eSBarry Smith /*
31c5083d92SDave May  default setup mode
321e07b27eSBarry Smith 
33c5083d92SDave May  [1a] scatter to (FORWARD)
341e07b27eSBarry Smith  x(comm) -> xtmp(comm)
35c5083d92SDave May  [1b] local copy (to) ranks with color = 0
361e07b27eSBarry Smith  xred(subcomm) <- xtmp
371e07b27eSBarry Smith 
38c5083d92SDave May  [2] solve on sub KSP to obtain yred(subcomm)
39c5083d92SDave May 
40c5083d92SDave May  [3a] local copy (from) ranks with color = 0
411e07b27eSBarry Smith  yred(subcomm) --> xtmp
42c5083d92SDave May  [2b] scatter from (REVERSE)
431e07b27eSBarry Smith  xtmp(comm) -> y(comm)
441e07b27eSBarry Smith */
451e07b27eSBarry Smith 
468d9f7141SDave May /*
478d9f7141SDave May   Collective on MPI_Comm[comm_f]
488d9f7141SDave May   Notes
498d9f7141SDave May    * Using comm_f = MPI_COMM_NULL will result in an error
508d9f7141SDave May    * Using comm_c = MPI_COMM_NULL is valid. If all instances of comm_c are NULL the subcomm is not valid.
518d9f7141SDave May    * If any non NULL comm_c communicator cannot map any of its ranks to comm_f, the subcomm is not valid.
528d9f7141SDave May */
538d9f7141SDave May PetscErrorCode PCTelescopeTestValidSubcomm(MPI_Comm comm_f,MPI_Comm comm_c,PetscBool *isvalid)
548d9f7141SDave May {
5557f12427SDave May   PetscInt       valid = 1;
568d9f7141SDave May   MPI_Group      group_f,group_c;
578d9f7141SDave May   PetscErrorCode ierr;
588d9f7141SDave May   PetscMPIInt    count,k,size_f = 0,size_c = 0,size_c_sum = 0;
598d9f7141SDave May   int            *ranks_f = NULL,*ranks_c = NULL;
608d9f7141SDave May 
6157f12427SDave May   PetscFunctionBegin;
628d9f7141SDave May   if (comm_f == MPI_COMM_NULL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"comm_f cannot be MPI_COMM_NULL");
638d9f7141SDave May 
648d9f7141SDave May   ierr = MPI_Comm_group(comm_f,&group_f);CHKERRQ(ierr);
658d9f7141SDave May   if (comm_c != MPI_COMM_NULL) {
668d9f7141SDave May     ierr = MPI_Comm_group(comm_c,&group_c);CHKERRQ(ierr);
678d9f7141SDave May   }
688d9f7141SDave May 
698d9f7141SDave May   ierr = MPI_Comm_size(comm_f,&size_f);CHKERRQ(ierr);
708d9f7141SDave May   if (comm_c != MPI_COMM_NULL) {
718d9f7141SDave May     ierr = MPI_Comm_size(comm_c,&size_c);CHKERRQ(ierr);
728d9f7141SDave May   }
738d9f7141SDave May 
748d9f7141SDave May   /* check not all comm_c's are NULL */
758d9f7141SDave May   size_c_sum = size_c;
768d9f7141SDave May   ierr = MPI_Allreduce(MPI_IN_PLACE,&size_c_sum,1,MPI_INT,MPI_SUM,comm_f);CHKERRQ(ierr);
778d9f7141SDave May   if (size_c_sum == 0) {
788d9f7141SDave May     valid = 0;
798d9f7141SDave May   }
808d9f7141SDave May 
818d9f7141SDave May   /* check we can map at least 1 rank in comm_c to comm_f */
828d9f7141SDave May   ierr = PetscMalloc1(size_f,&ranks_f);CHKERRQ(ierr);
838d9f7141SDave May   ierr = PetscMalloc1(size_c,&ranks_c);CHKERRQ(ierr);
848d9f7141SDave May   for (k=0; k<size_f; k++) {
858d9f7141SDave May     ranks_f[k] = MPI_UNDEFINED;
868d9f7141SDave May   }
878d9f7141SDave May   for (k=0; k<size_c; k++) {
888d9f7141SDave May     ranks_c[k] = (int)k;
898d9f7141SDave May   }
908d9f7141SDave May 
9157f12427SDave May   /*
9257f12427SDave May    MPI_Group_translate_ranks() returns a non-zero exit code if any rank cannot be translated.
9357f12427SDave May    I do not want the code to terminate immediately if this occurs, rather I want to throw
9457f12427SDave May    the error later (during PCSetUp_Telescope()) via SETERRQ() with a message indicating
9557f12427SDave May    that comm_c is not a valid sub-communicator.
9657f12427SDave May    Hence I purposefully do not call CHKERRQ() after MPI_Group_translate_ranks().
9757f12427SDave May   */
988d9f7141SDave May   count = 0;
998d9f7141SDave May   if (comm_c != MPI_COMM_NULL) {
10066b79024SDave May     (void)MPI_Group_translate_ranks(group_c,size_c,ranks_c,group_f,ranks_f);
1018d9f7141SDave May     for (k=0; k<size_f; k++) {
1028d9f7141SDave May       if (ranks_f[k] == MPI_UNDEFINED) {
1038d9f7141SDave May         count++;
1048d9f7141SDave May       }
1058d9f7141SDave May     }
1068d9f7141SDave May   }
1078d9f7141SDave May   if (count == size_f) {
1088d9f7141SDave May     valid = 0;
1098d9f7141SDave May   }
1108d9f7141SDave May 
11157f12427SDave May   ierr = MPI_Allreduce(MPI_IN_PLACE,&valid,1,MPIU_INT,MPI_MIN,comm_f);CHKERRQ(ierr);
1128d9f7141SDave May   if (valid == 1) { *isvalid = PETSC_TRUE; }
1138d9f7141SDave May   else { *isvalid = PETSC_FALSE; }
1148d9f7141SDave May 
1158d9f7141SDave May   ierr = PetscFree(ranks_f);CHKERRQ(ierr);
1168d9f7141SDave May   ierr = PetscFree(ranks_c);CHKERRQ(ierr);
1178d9f7141SDave May   ierr = MPI_Group_free(&group_f);CHKERRQ(ierr);
1188d9f7141SDave May   if (comm_c != MPI_COMM_NULL) {
1198d9f7141SDave May     ierr = MPI_Group_free(&group_c);CHKERRQ(ierr);
1208d9f7141SDave May   }
1218d9f7141SDave May   PetscFunctionReturn(0);
1228d9f7141SDave May }
1238d9f7141SDave May 
1241e07b27eSBarry Smith DM private_PCTelescopeGetSubDM(PC_Telescope sred)
1251e07b27eSBarry Smith {
126c6a0d831SBarry Smith   DM subdm = NULL;
1271e07b27eSBarry Smith 
12857f12427SDave May   if (!PCTelescope_isActiveRank(sred)) { subdm = NULL; }
1291e07b27eSBarry Smith   else {
1301e07b27eSBarry Smith     switch (sred->sr_type) {
1311e07b27eSBarry Smith     case TELESCOPE_DEFAULT: subdm = NULL;
1321e07b27eSBarry Smith       break;
1331e07b27eSBarry Smith     case TELESCOPE_DMDA:    subdm = ((PC_Telescope_DMDACtx*)sred->dm_ctx)->dmrepart;
1341e07b27eSBarry Smith       break;
1351e07b27eSBarry Smith     case TELESCOPE_DMPLEX:  subdm = NULL;
1361e07b27eSBarry Smith       break;
1378d9f7141SDave May     case TELESCOPE_COARSEDM: if (sred->ksp) { KSPGetDM(sred->ksp,&subdm); }
1388d9f7141SDave May       break;
1391e07b27eSBarry Smith     }
1401e07b27eSBarry Smith   }
1411e07b27eSBarry Smith   return(subdm);
1421e07b27eSBarry Smith }
1431e07b27eSBarry Smith 
1441e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_default(PC pc,PC_Telescope sred)
1451e07b27eSBarry Smith {
1461e07b27eSBarry Smith   PetscErrorCode ierr;
1471e07b27eSBarry Smith   PetscInt       m,M,bs,st,ed;
1481e07b27eSBarry Smith   Vec            x,xred,yred,xtmp;
1491e07b27eSBarry Smith   Mat            B;
1501e07b27eSBarry Smith   MPI_Comm       comm,subcomm;
1511e07b27eSBarry Smith   VecScatter     scatter;
1521e07b27eSBarry Smith   IS             isin;
1531e07b27eSBarry Smith 
1541e07b27eSBarry Smith   PetscFunctionBegin;
1551e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: setup (default)\n");CHKERRQ(ierr);
1561e07b27eSBarry Smith   comm = PetscSubcommParent(sred->psubcomm);
1571e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
1581e07b27eSBarry Smith 
1591e07b27eSBarry Smith   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
1601e07b27eSBarry Smith   ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr);
1611e07b27eSBarry Smith   ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr);
1621e07b27eSBarry Smith   ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr);
1631e07b27eSBarry Smith 
1641e07b27eSBarry Smith   xred = NULL;
1653ac26c5eSBarry Smith   m    = 0;
16657f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
1671e07b27eSBarry Smith     ierr = VecCreate(subcomm,&xred);CHKERRQ(ierr);
1681e07b27eSBarry Smith     ierr = VecSetSizes(xred,PETSC_DECIDE,M);CHKERRQ(ierr);
1691e07b27eSBarry Smith     ierr = VecSetBlockSize(xred,bs);CHKERRQ(ierr);
1701e07b27eSBarry Smith     ierr = VecSetFromOptions(xred);CHKERRQ(ierr);
171ca43db0aSBarry Smith     ierr = VecGetLocalSize(xred,&m);CHKERRQ(ierr);
1721e07b27eSBarry Smith   }
1731e07b27eSBarry Smith 
1741e07b27eSBarry Smith   yred = NULL;
17557f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
1761e07b27eSBarry Smith     ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr);
1771e07b27eSBarry Smith   }
1781e07b27eSBarry Smith 
1791e07b27eSBarry Smith   ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr);
1801e07b27eSBarry Smith   ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr);
1811e07b27eSBarry Smith   ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr);
1821e07b27eSBarry Smith   ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr);
1831e07b27eSBarry Smith 
18457f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
1851e07b27eSBarry Smith     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
1861e07b27eSBarry Smith     ierr = ISCreateStride(comm,(ed-st),st,1,&isin);CHKERRQ(ierr);
1871e07b27eSBarry Smith   } else {
1881e07b27eSBarry Smith     ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr);
1893ac26c5eSBarry Smith     ierr = ISCreateStride(comm,0,st,1,&isin);CHKERRQ(ierr);
1901e07b27eSBarry Smith   }
1911e07b27eSBarry Smith   ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr);
1921e07b27eSBarry Smith 
1939448b7f1SJunchao Zhang   ierr = VecScatterCreate(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr);
1941e07b27eSBarry Smith 
1951e07b27eSBarry Smith   sred->isin    = isin;
1961e07b27eSBarry Smith   sred->scatter = scatter;
1971e07b27eSBarry Smith   sred->xred    = xred;
1981e07b27eSBarry Smith   sred->yred    = yred;
1991e07b27eSBarry Smith   sred->xtmp    = xtmp;
2001e07b27eSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
2011e07b27eSBarry Smith   PetscFunctionReturn(0);
2021e07b27eSBarry Smith }
2031e07b27eSBarry Smith 
2041e07b27eSBarry Smith PetscErrorCode PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A)
2051e07b27eSBarry Smith {
2061e07b27eSBarry Smith   PetscErrorCode ierr;
2071e07b27eSBarry Smith   MPI_Comm       comm,subcomm;
2081e07b27eSBarry Smith   Mat            Bred,B;
2091e07b27eSBarry Smith   PetscInt       nr,nc;
2101e07b27eSBarry Smith   IS             isrow,iscol;
2111e07b27eSBarry Smith   Mat            Blocal,*_Blocal;
2121e07b27eSBarry Smith 
2131e07b27eSBarry Smith   PetscFunctionBegin;
2141e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (default)\n");CHKERRQ(ierr);
2151e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
2161e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
2171e07b27eSBarry Smith   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
2181e07b27eSBarry Smith   ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr);
2191e07b27eSBarry Smith   isrow = sred->isin;
2201e07b27eSBarry Smith   ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr);
2217dae84e0SHong Zhang   ierr = MatCreateSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr);
2221e07b27eSBarry Smith   Blocal = *_Blocal;
2231e07b27eSBarry Smith   ierr = PetscFree(_Blocal);CHKERRQ(ierr);
2241e07b27eSBarry Smith   Bred = NULL;
22557f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
2261e07b27eSBarry Smith     PetscInt mm;
2271e07b27eSBarry Smith 
2281e07b27eSBarry Smith     if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; }
2291e07b27eSBarry Smith 
2301e07b27eSBarry Smith     ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr);
2311e07b27eSBarry Smith     ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr);
2321e07b27eSBarry Smith   }
2331e07b27eSBarry Smith   *A = Bred;
2341e07b27eSBarry Smith   ierr = ISDestroy(&iscol);CHKERRQ(ierr);
2351e07b27eSBarry Smith   ierr = MatDestroy(&Blocal);CHKERRQ(ierr);
2361e07b27eSBarry Smith   PetscFunctionReturn(0);
2371e07b27eSBarry Smith }
2381e07b27eSBarry Smith 
239392968a1SPatrick Sanan static PetscErrorCode PCTelescopeSubNullSpaceCreate_Telescope(PC pc,PC_Telescope sred,MatNullSpace nullspace,MatNullSpace *sub_nullspace)
2401e07b27eSBarry Smith {
2411e07b27eSBarry Smith   PetscErrorCode ierr;
2421e07b27eSBarry Smith   PetscBool      has_const;
2431e07b27eSBarry Smith   const Vec      *vecs;
244c41e779fSDave May   Vec            *sub_vecs = NULL;
245392968a1SPatrick Sanan   PetscInt       i,k,n = 0;
2461e07b27eSBarry Smith   MPI_Comm       subcomm;
2471e07b27eSBarry Smith 
2481e07b27eSBarry Smith   PetscFunctionBegin;
2491e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
2501e07b27eSBarry Smith   ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr);
2511e07b27eSBarry Smith 
25257f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
253e3acf2f7SBarry Smith     if (n) {
254e3acf2f7SBarry Smith       ierr = VecDuplicateVecs(sred->xred,n,&sub_vecs);CHKERRQ(ierr);
2551e07b27eSBarry Smith     }
2561e07b27eSBarry Smith   }
2571e07b27eSBarry Smith 
2581e07b27eSBarry Smith   /* copy entries */
2591e07b27eSBarry Smith   for (k=0; k<n; k++) {
2601e07b27eSBarry Smith     const PetscScalar *x_array;
2611e07b27eSBarry Smith     PetscScalar       *LA_sub_vec;
26213c30530SDave May     PetscInt          st,ed;
2631e07b27eSBarry Smith 
2641e07b27eSBarry Smith     /* pull in vector x->xtmp */
2651e07b27eSBarry Smith     ierr = VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2661e07b27eSBarry Smith     ierr = VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
26747856c66SBarry Smith     if (sub_vecs) {
268a04a6428SPatrick Sanan       /* copy vector entries into xred */
2691e07b27eSBarry Smith       ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
270ea2b237eSDave May       if (sub_vecs[k]) {
2711e07b27eSBarry Smith         ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr);
2721e07b27eSBarry Smith         ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
2731e07b27eSBarry Smith         for (i=0; i<ed-st; i++) {
2741e07b27eSBarry Smith           LA_sub_vec[i] = x_array[i];
2751e07b27eSBarry Smith         }
2761e07b27eSBarry Smith         ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
2771e07b27eSBarry Smith       }
2781e07b27eSBarry Smith       ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
2791e07b27eSBarry Smith     }
28047856c66SBarry Smith   }
2811e07b27eSBarry Smith 
28257f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
283d8b9d5b7SPatrick Sanan     /* create new (near) nullspace for redundant object */
284392968a1SPatrick Sanan     ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,sub_nullspace);CHKERRQ(ierr);
285392968a1SPatrick Sanan     ierr = VecDestroyVecs(n,&sub_vecs);CHKERRQ(ierr);
286d8b9d5b7SPatrick Sanan     if (nullspace->remove) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope");
287d8b9d5b7SPatrick Sanan     if (nullspace->rmctx) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callback context not supported when propagating (near) nullspaces with PCTelescope");
288d8b9d5b7SPatrick Sanan   }
289392968a1SPatrick Sanan   PetscFunctionReturn(0);
290392968a1SPatrick Sanan }
291392968a1SPatrick Sanan 
292392968a1SPatrick Sanan static PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat)
293392968a1SPatrick Sanan {
294392968a1SPatrick Sanan   PetscErrorCode ierr;
295392968a1SPatrick Sanan   Mat            B;
296392968a1SPatrick Sanan 
297392968a1SPatrick Sanan   PetscFunctionBegin;
298392968a1SPatrick Sanan   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
299392968a1SPatrick Sanan   /* Propagate the nullspace if it exists */
300392968a1SPatrick Sanan   {
301392968a1SPatrick Sanan     MatNullSpace nullspace,sub_nullspace;
302392968a1SPatrick Sanan     ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr);
303392968a1SPatrick Sanan     if (nullspace) {
304392968a1SPatrick Sanan       ierr = PetscInfo(pc,"PCTelescope: generating nullspace (default)\n");CHKERRQ(ierr);
305392968a1SPatrick Sanan       ierr = PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nullspace,&sub_nullspace);CHKERRQ(ierr);
30657f12427SDave May       if (PCTelescope_isActiveRank(sred)) {
307392968a1SPatrick Sanan         ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr);
30841ff1ee9SPatrick Sanan         ierr = MatNullSpaceDestroy(&sub_nullspace);CHKERRQ(ierr);
3091e07b27eSBarry Smith       }
310392968a1SPatrick Sanan     }
311392968a1SPatrick Sanan   }
312392968a1SPatrick Sanan   /* Propagate the near nullspace if it exists */
313392968a1SPatrick Sanan   {
314392968a1SPatrick Sanan     MatNullSpace nearnullspace,sub_nearnullspace;
315392968a1SPatrick Sanan     ierr = MatGetNearNullSpace(B,&nearnullspace);CHKERRQ(ierr);
316392968a1SPatrick Sanan     if (nearnullspace) {
317392968a1SPatrick Sanan       ierr = PetscInfo(pc,"PCTelescope: generating near nullspace (default)\n");CHKERRQ(ierr);
318392968a1SPatrick Sanan       ierr = PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nearnullspace,&sub_nearnullspace);CHKERRQ(ierr);
31957f12427SDave May       if (PCTelescope_isActiveRank(sred)) {
320392968a1SPatrick Sanan         ierr = MatSetNearNullSpace(sub_mat,sub_nearnullspace);CHKERRQ(ierr);
321392968a1SPatrick Sanan         ierr = MatNullSpaceDestroy(&sub_nearnullspace);CHKERRQ(ierr);
322392968a1SPatrick Sanan       }
323392968a1SPatrick Sanan     }
324392968a1SPatrick Sanan   }
3251e07b27eSBarry Smith   PetscFunctionReturn(0);
3261e07b27eSBarry Smith }
3271e07b27eSBarry Smith 
3281e07b27eSBarry Smith static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer)
3291e07b27eSBarry Smith {
3301e07b27eSBarry Smith   PC_Telescope   sred = (PC_Telescope)pc->data;
3311e07b27eSBarry Smith   PetscErrorCode ierr;
3321e07b27eSBarry Smith   PetscBool      iascii,isstring;
3331e07b27eSBarry Smith   PetscViewer    subviewer;
3341e07b27eSBarry Smith 
3351e07b27eSBarry Smith   PetscFunctionBegin;
3361e07b27eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
3371e07b27eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
3381e07b27eSBarry Smith   if (iascii) {
3398d9f7141SDave May     {
3401e07b27eSBarry Smith       MPI_Comm    comm,subcomm;
3411e07b27eSBarry Smith       PetscMPIInt comm_size,subcomm_size;
3428d9f7141SDave May       DM          dm = NULL,subdm = NULL;
3431e07b27eSBarry Smith 
3441e07b27eSBarry Smith       ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
3451e07b27eSBarry Smith       subdm = private_PCTelescopeGetSubDM(sred);
3468d9f7141SDave May 
3478d9f7141SDave May       if (sred->psubcomm) {
3481e07b27eSBarry Smith         comm = PetscSubcommParent(sred->psubcomm);
3491e07b27eSBarry Smith         subcomm = PetscSubcommChild(sred->psubcomm);
3501e07b27eSBarry Smith         ierr = MPI_Comm_size(comm,&comm_size);CHKERRQ(ierr);
3511e07b27eSBarry Smith         ierr = MPI_Comm_size(subcomm,&subcomm_size);CHKERRQ(ierr);
3521e07b27eSBarry Smith 
3538d9f7141SDave May         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3548d9f7141SDave May         ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm: parent comm size reduction factor = %D\n",sred->redfactor);CHKERRQ(ierr);
3558d9f7141SDave May         ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm: parent_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size);CHKERRQ(ierr);
35648a10b22SPatrick Sanan         switch (sred->subcommtype) {
35748a10b22SPatrick Sanan         case PETSC_SUBCOMM_INTERLACED :
3588d9f7141SDave May           ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm: type = interlaced\n",sred->subcommtype);CHKERRQ(ierr);
35948a10b22SPatrick Sanan           break;
36048a10b22SPatrick Sanan         case PETSC_SUBCOMM_CONTIGUOUS :
3618d9f7141SDave May           ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm type = contiguous\n",sred->subcommtype);CHKERRQ(ierr);
36248a10b22SPatrick Sanan           break;
36348a10b22SPatrick Sanan         default :
36448a10b22SPatrick Sanan           SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"General subcomm type not supported by PCTelescope");
36548a10b22SPatrick Sanan         }
3668d9f7141SDave May         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3678d9f7141SDave May       } else {
3688d9f7141SDave May         ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
3698d9f7141SDave May         subcomm = sred->subcomm;
37057f12427SDave May         if (!PCTelescope_isActiveRank(sred)) {
3718d9f7141SDave May           subcomm = PETSC_COMM_SELF;
3728d9f7141SDave May         }
3738d9f7141SDave May 
3748d9f7141SDave May         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3758d9f7141SDave May         ierr = PetscViewerASCIIPrintf(viewer,"subcomm: using user provided sub-communicator\n");CHKERRQ(ierr);
3768d9f7141SDave May         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3778d9f7141SDave May       }
3788d9f7141SDave May 
3791e07b27eSBarry Smith       ierr = PetscViewerGetSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
38057f12427SDave May       if (PCTelescope_isActiveRank(sred)) {
3811e07b27eSBarry Smith         ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr);
3821e07b27eSBarry Smith 
3831e07b27eSBarry Smith         if (dm && sred->ignore_dm) {
384efd4aadfSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"ignoring DM\n");CHKERRQ(ierr);
3851e07b27eSBarry Smith         }
3867c5279cbSDave May         if (sred->ignore_kspcomputeoperators) {
387efd4aadfSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"ignoring KSPComputeOperators\n");CHKERRQ(ierr);
3887c5279cbSDave May         }
3891e07b27eSBarry Smith         switch (sred->sr_type) {
3901e07b27eSBarry Smith         case TELESCOPE_DEFAULT:
3918d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"setup type: default\n");CHKERRQ(ierr);
3921e07b27eSBarry Smith           break;
3931e07b27eSBarry Smith         case TELESCOPE_DMDA:
3948d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"setup type: DMDA auto-repartitioning\n");CHKERRQ(ierr);
3958ef9ca65SPatrick Sanan           ierr = DMView_DA_Short(subdm,subviewer);CHKERRQ(ierr);
3961e07b27eSBarry Smith           break;
3971e07b27eSBarry Smith         case TELESCOPE_DMPLEX:
3988d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"setup type: DMPLEX auto-repartitioning\n");CHKERRQ(ierr);
3991e07b27eSBarry Smith           break;
4008d9f7141SDave May         case TELESCOPE_COARSEDM:
4018d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"setup type: coarse DM\n");CHKERRQ(ierr);
4028d9f7141SDave May           break;
4038d9f7141SDave May         }
4048d9f7141SDave May 
4058d9f7141SDave May         if (dm) {
4068d9f7141SDave May           PetscObject obj = (PetscObject)dm;
4078d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"Parent DM object:");CHKERRQ(ierr);
4088d9f7141SDave May           PetscViewerASCIIUseTabs(subviewer,PETSC_FALSE);
4098d9f7141SDave May           if (obj->type_name) { PetscViewerASCIIPrintf(subviewer," type = %s;",obj->type_name); }
4108d9f7141SDave May           if (obj->name) { PetscViewerASCIIPrintf(subviewer," name = %s;",obj->name); }
4118d9f7141SDave May           if (obj->prefix) { PetscViewerASCIIPrintf(subviewer," prefix = %s",obj->prefix); }
4128d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"\n");CHKERRQ(ierr);
4138d9f7141SDave May           PetscViewerASCIIUseTabs(subviewer,PETSC_TRUE);
4148d9f7141SDave May         } else {
4158d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"Parent DM object: NULL\n");CHKERRQ(ierr);
4168d9f7141SDave May         }
4178d9f7141SDave May         if (subdm) {
4188d9f7141SDave May           PetscObject obj = (PetscObject)subdm;
4198d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"Sub DM object:");CHKERRQ(ierr);
4208d9f7141SDave May           PetscViewerASCIIUseTabs(subviewer,PETSC_FALSE);
4218d9f7141SDave May           if (obj->type_name) { PetscViewerASCIIPrintf(subviewer," type = %s;",obj->type_name); }
4228d9f7141SDave May           if (obj->name) { PetscViewerASCIIPrintf(subviewer," name = %s;",obj->name); }
4238d9f7141SDave May           if (obj->prefix) { PetscViewerASCIIPrintf(subviewer," prefix = %s",obj->prefix); }
4248d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"\n");CHKERRQ(ierr);
4258d9f7141SDave May           PetscViewerASCIIUseTabs(subviewer,PETSC_TRUE);
4268d9f7141SDave May         } else {
4278d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"Sub DM object: NULL\n");CHKERRQ(ierr);
4281e07b27eSBarry Smith         }
4291e07b27eSBarry Smith 
4301e07b27eSBarry Smith         ierr = KSPView(sred->ksp,subviewer);CHKERRQ(ierr);
4311e07b27eSBarry Smith         ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr);
4321e07b27eSBarry Smith       }
4331e07b27eSBarry Smith       ierr = PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
4341e07b27eSBarry Smith     }
4351e07b27eSBarry Smith   }
4361e07b27eSBarry Smith   PetscFunctionReturn(0);
4371e07b27eSBarry Smith }
4381e07b27eSBarry Smith 
4391e07b27eSBarry Smith static PetscErrorCode PCSetUp_Telescope(PC pc)
4401e07b27eSBarry Smith {
4411e07b27eSBarry Smith   PC_Telescope    sred = (PC_Telescope)pc->data;
4421e07b27eSBarry Smith   PetscErrorCode  ierr;
443bd49479cSSatish Balay   MPI_Comm        comm,subcomm=0;
4441e07b27eSBarry Smith   PCTelescopeType sr_type;
4451e07b27eSBarry Smith 
4461e07b27eSBarry Smith   PetscFunctionBegin;
4471e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4481e07b27eSBarry Smith 
4491e07b27eSBarry Smith   /* Determine type of setup/update */
4501e07b27eSBarry Smith   if (!pc->setupcalled) {
4511e07b27eSBarry Smith     PetscBool has_dm,same;
4521e07b27eSBarry Smith     DM        dm;
4531e07b27eSBarry Smith 
4541e07b27eSBarry Smith     sr_type = TELESCOPE_DEFAULT;
4551e07b27eSBarry Smith     has_dm = PETSC_FALSE;
4561e07b27eSBarry Smith     ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
4571e07b27eSBarry Smith     if (dm) { has_dm = PETSC_TRUE; }
4581e07b27eSBarry Smith     if (has_dm) {
4591e07b27eSBarry Smith       /* check for dmda */
4601e07b27eSBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&same);CHKERRQ(ierr);
4611e07b27eSBarry Smith       if (same) {
4621e07b27eSBarry Smith         ierr = PetscInfo(pc,"PCTelescope: found DMDA\n");CHKERRQ(ierr);
4631e07b27eSBarry Smith         sr_type = TELESCOPE_DMDA;
4641e07b27eSBarry Smith       }
4651e07b27eSBarry Smith       /* check for dmplex */
4661e07b27eSBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same);CHKERRQ(ierr);
4671e07b27eSBarry Smith       if (same) {
468994fe344SLisandro Dalcin         ierr = PetscInfo(pc,"PCTelescope: found DMPLEX\n");CHKERRQ(ierr);
4691e07b27eSBarry Smith         sr_type = TELESCOPE_DMPLEX;
4701e07b27eSBarry Smith       }
4718d9f7141SDave May 
4728d9f7141SDave May       if (sred->use_coarse_dm) {
4738d9f7141SDave May         ierr = PetscInfo(pc,"PCTelescope: using coarse DM\n");CHKERRQ(ierr);
4748d9f7141SDave May         sr_type = TELESCOPE_COARSEDM;
4751e07b27eSBarry Smith       }
4761e07b27eSBarry Smith 
4771e07b27eSBarry Smith       if (sred->ignore_dm) {
4788d9f7141SDave May         ierr = PetscInfo(pc,"PCTelescope: ignoring DM\n");CHKERRQ(ierr);
4791e07b27eSBarry Smith         sr_type = TELESCOPE_DEFAULT;
4801e07b27eSBarry Smith       }
4818d9f7141SDave May     }
4821e07b27eSBarry Smith     sred->sr_type = sr_type;
4831e07b27eSBarry Smith   } else {
4841e07b27eSBarry Smith     sr_type = sred->sr_type;
4851e07b27eSBarry Smith   }
4861e07b27eSBarry Smith 
487d8b9d5b7SPatrick Sanan   /* set function pointers for repartition setup, matrix creation/update, matrix (near) nullspace, and reset functionality */
4881e07b27eSBarry Smith   switch (sr_type) {
4891e07b27eSBarry Smith   case TELESCOPE_DEFAULT:
4901e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
4911e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
4921e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
4931e07b27eSBarry Smith     sred->pctelescope_reset_type              = NULL;
4941e07b27eSBarry Smith     break;
4951e07b27eSBarry Smith   case TELESCOPE_DMDA:
4961e07b27eSBarry Smith     pc->ops->apply                            = PCApply_Telescope_dmda;
497f650675bSDave May     pc->ops->applyrichardson                  = PCApplyRichardson_Telescope_dmda;
4981e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_dmda;
4991e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_dmda;
5001e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda;
5011e07b27eSBarry Smith     sred->pctelescope_reset_type              = PCReset_Telescope_dmda;
5021e07b27eSBarry Smith     break;
503d8b9d5b7SPatrick Sanan   case TELESCOPE_DMPLEX: SETERRQ(comm,PETSC_ERR_SUP,"Support for DMPLEX is currently not available");
5041e07b27eSBarry Smith     break;
5058d9f7141SDave May   case TELESCOPE_COARSEDM:
5068d9f7141SDave May     pc->ops->apply                            = PCApply_Telescope_CoarseDM;
5078d9f7141SDave May     pc->ops->applyrichardson                  = PCApplyRichardson_Telescope_CoarseDM;
5088d9f7141SDave May     sred->pctelescope_setup_type              = PCTelescopeSetUp_CoarseDM;
5098d9f7141SDave May     sred->pctelescope_matcreate_type          = NULL;
5108d9f7141SDave May     sred->pctelescope_matnullspacecreate_type = NULL; /* PCTelescopeMatNullSpaceCreate_CoarseDM; */
5118d9f7141SDave May     sred->pctelescope_reset_type              = PCReset_Telescope_CoarseDM;
5121e07b27eSBarry Smith     break;
5138d9f7141SDave May   default: SETERRQ(comm,PETSC_ERR_SUP,"Support only provided for: repartitioning an operator; repartitioning a DMDA; or using a coarse DM");
5148d9f7141SDave May     break;
5158d9f7141SDave May   }
5168d9f7141SDave May 
5178d9f7141SDave May   /* subcomm definition */
5188d9f7141SDave May   if (!pc->setupcalled) {
5198d9f7141SDave May     if ((sr_type == TELESCOPE_DEFAULT) || (sr_type == TELESCOPE_DMDA)) {
5208d9f7141SDave May       if (!sred->psubcomm) {
5218d9f7141SDave May         ierr = PetscSubcommCreate(comm,&sred->psubcomm);CHKERRQ(ierr);
5228d9f7141SDave May         ierr = PetscSubcommSetNumber(sred->psubcomm,sred->redfactor);CHKERRQ(ierr);
5238d9f7141SDave May         ierr = PetscSubcommSetType(sred->psubcomm,sred->subcommtype);CHKERRQ(ierr);
5248d9f7141SDave May         ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr);
5258d9f7141SDave May         sred->subcomm = PetscSubcommChild(sred->psubcomm);
5268d9f7141SDave May       }
5278d9f7141SDave May     } else { /* query PC for DM, check communicators */
5288d9f7141SDave May       DM          dm,dm_coarse_partition = NULL;
5298d9f7141SDave May       MPI_Comm    comm_fine,comm_coarse_partition = MPI_COMM_NULL;
5308d9f7141SDave May       PetscMPIInt csize_fine=0,csize_coarse_partition=0,cs[2],csg[2],cnt=0;
5318d9f7141SDave May       PetscBool   isvalidsubcomm;
5328d9f7141SDave May 
5338d9f7141SDave May       ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
5348d9f7141SDave May       comm_fine = PetscObjectComm((PetscObject)dm);
5358d9f7141SDave May       ierr = DMGetCoarseDM(dm,&dm_coarse_partition);CHKERRQ(ierr);
5368d9f7141SDave May       if (dm_coarse_partition) { cnt = 1; }
5378d9f7141SDave May       ierr = MPI_Allreduce(MPI_IN_PLACE,&cnt,1,MPI_INT,MPI_SUM,comm_fine);CHKERRQ(ierr);
5388d9f7141SDave May       if (cnt == 0) SETERRQ(comm_fine,PETSC_ERR_SUP,"Zero instances of a coarse DM were found");
5398d9f7141SDave May 
5408d9f7141SDave May       ierr = MPI_Comm_size(comm_fine,&csize_fine);CHKERRQ(ierr);
5418d9f7141SDave May       if (dm_coarse_partition) {
5428d9f7141SDave May         comm_coarse_partition = PetscObjectComm((PetscObject)dm_coarse_partition);
5438d9f7141SDave May         ierr = MPI_Comm_size(comm_coarse_partition,&csize_coarse_partition);CHKERRQ(ierr);
5448d9f7141SDave May       }
5458d9f7141SDave May 
5468d9f7141SDave May       cs[0] = csize_fine;
5478d9f7141SDave May       cs[1] = csize_coarse_partition;
5488d9f7141SDave May       ierr = MPI_Allreduce(cs,csg,2,MPI_INT,MPI_MAX,comm_fine);CHKERRQ(ierr);
5498d9f7141SDave May       if (csg[0] == csg[1]) SETERRQ(comm_fine,PETSC_ERR_SUP,"Coarse DM uses the same size communicator as the parent DM attached to the PC");
5508d9f7141SDave May 
5518d9f7141SDave May       ierr = PCTelescopeTestValidSubcomm(comm_fine,comm_coarse_partition,&isvalidsubcomm);CHKERRQ(ierr);
5528d9f7141SDave May       if (!isvalidsubcomm) SETERRQ(comm_fine,PETSC_ERR_SUP,"Coarse DM communicator is not a sub-communicator of parentDM->comm");
5538d9f7141SDave May       sred->subcomm = comm_coarse_partition;
5548d9f7141SDave May     }
5558d9f7141SDave May   }
5568d9f7141SDave May   subcomm = sred->subcomm;
5578d9f7141SDave May 
5588d9f7141SDave May   /* internal KSP */
5598d9f7141SDave May   if (!pc->setupcalled) {
5608d9f7141SDave May     const char *prefix;
5618d9f7141SDave May 
56257f12427SDave May     if (PCTelescope_isActiveRank(sred)) {
5638d9f7141SDave May       ierr = KSPCreate(subcomm,&sred->ksp);CHKERRQ(ierr);
5648d9f7141SDave May       ierr = KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure);CHKERRQ(ierr);
5658d9f7141SDave May       ierr = PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
5668d9f7141SDave May       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp);CHKERRQ(ierr);
5678d9f7141SDave May       ierr = KSPSetType(sred->ksp,KSPPREONLY);CHKERRQ(ierr);
5688d9f7141SDave May       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
5698d9f7141SDave May       ierr = KSPSetOptionsPrefix(sred->ksp,prefix);CHKERRQ(ierr);
5708d9f7141SDave May       ierr = KSPAppendOptionsPrefix(sred->ksp,"telescope_");CHKERRQ(ierr);
5718d9f7141SDave May     }
5721e07b27eSBarry Smith   }
5731e07b27eSBarry Smith 
5741e07b27eSBarry Smith   /* setup */
575aaa7a805SDave May   if (!pc->setupcalled && sred->pctelescope_setup_type) {
5761e07b27eSBarry Smith     ierr = sred->pctelescope_setup_type(pc,sred);CHKERRQ(ierr);
5771e07b27eSBarry Smith   }
5781e07b27eSBarry Smith   /* update */
5791e07b27eSBarry Smith   if (!pc->setupcalled) {
5801e07b27eSBarry Smith     if (sred->pctelescope_matcreate_type) {
5811e07b27eSBarry Smith       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred);CHKERRQ(ierr);
5821e07b27eSBarry Smith     }
5831e07b27eSBarry Smith     if (sred->pctelescope_matnullspacecreate_type) {
584392968a1SPatrick Sanan       ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr);
5851e07b27eSBarry Smith     }
5861e07b27eSBarry Smith   } else {
5871e07b27eSBarry Smith     if (sred->pctelescope_matcreate_type) {
5881e07b27eSBarry Smith       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred);CHKERRQ(ierr);
5891e07b27eSBarry Smith     }
5901e07b27eSBarry Smith   }
5911e07b27eSBarry Smith 
5921e07b27eSBarry Smith   /* common - no construction */
59357f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
5941e07b27eSBarry Smith     ierr = KSPSetOperators(sred->ksp,sred->Bred,sred->Bred);CHKERRQ(ierr);
5951e07b27eSBarry Smith     if (pc->setfromoptionscalled && !pc->setupcalled){
5961e07b27eSBarry Smith       ierr = KSPSetFromOptions(sred->ksp);CHKERRQ(ierr);
5971e07b27eSBarry Smith     }
5981e07b27eSBarry Smith   }
5991e07b27eSBarry Smith   PetscFunctionReturn(0);
6001e07b27eSBarry Smith }
6011e07b27eSBarry Smith 
6021e07b27eSBarry Smith static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y)
6031e07b27eSBarry Smith {
6041e07b27eSBarry Smith   PC_Telescope      sred = (PC_Telescope)pc->data;
6051e07b27eSBarry Smith   PetscErrorCode    ierr;
6061e07b27eSBarry Smith   Vec               xtmp,xred,yred;
60713c30530SDave May   PetscInt          i,st,ed;
6081e07b27eSBarry Smith   VecScatter        scatter;
6091e07b27eSBarry Smith   PetscScalar       *array;
6101e07b27eSBarry Smith   const PetscScalar *x_array;
6111e07b27eSBarry Smith 
6121e07b27eSBarry Smith   PetscFunctionBegin;
613bf00f589SPatrick Sanan   ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr);
614bf00f589SPatrick Sanan 
6151e07b27eSBarry Smith   xtmp    = sred->xtmp;
6161e07b27eSBarry Smith   scatter = sred->scatter;
6171e07b27eSBarry Smith   xred    = sred->xred;
6181e07b27eSBarry Smith   yred    = sred->yred;
6191e07b27eSBarry Smith 
6201e07b27eSBarry Smith   /* pull in vector x->xtmp */
6211e07b27eSBarry Smith   ierr = VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6221e07b27eSBarry Smith   ierr = VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6231e07b27eSBarry Smith 
624bf00f589SPatrick Sanan   /* copy vector entries into xred */
6251e07b27eSBarry Smith   ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr);
6261e07b27eSBarry Smith   if (xred) {
6271e07b27eSBarry Smith     PetscScalar *LA_xred;
6281e07b27eSBarry Smith     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
6291e07b27eSBarry Smith     ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr);
6301e07b27eSBarry Smith     for (i=0; i<ed-st; i++) {
6311e07b27eSBarry Smith       LA_xred[i] = x_array[i];
6321e07b27eSBarry Smith     }
6331e07b27eSBarry Smith     ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr);
6341e07b27eSBarry Smith   }
6351e07b27eSBarry Smith   ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr);
6361e07b27eSBarry Smith   /* solve */
63757f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
6381e07b27eSBarry Smith     ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr);
639c0decd05SBarry Smith     ierr = KSPCheckSolve(sred->ksp,pc,yred);CHKERRQ(ierr);
6401e07b27eSBarry Smith   }
6411e07b27eSBarry Smith   /* return vector */
6421e07b27eSBarry Smith   ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr);
6431e07b27eSBarry Smith   if (yred) {
6441e07b27eSBarry Smith     const PetscScalar *LA_yred;
6451e07b27eSBarry Smith     ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr);
6461e07b27eSBarry Smith     ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr);
6471e07b27eSBarry Smith     for (i=0; i<ed-st; i++) {
6481e07b27eSBarry Smith       array[i] = LA_yred[i];
6491e07b27eSBarry Smith     }
6501e07b27eSBarry Smith     ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr);
6511e07b27eSBarry Smith   }
6521e07b27eSBarry Smith   ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr);
6531e07b27eSBarry Smith   ierr = VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6541e07b27eSBarry Smith   ierr = VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6551e07b27eSBarry Smith   PetscFunctionReturn(0);
6561e07b27eSBarry Smith }
6571e07b27eSBarry Smith 
658f650675bSDave May 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)
659f650675bSDave May {
660f650675bSDave May   PC_Telescope      sred = (PC_Telescope)pc->data;
661f650675bSDave May   PetscErrorCode    ierr;
662a1d91a28SDave May   Vec               xtmp,yred;
663f650675bSDave May   PetscInt          i,st,ed;
664f650675bSDave May   VecScatter        scatter;
665f650675bSDave May   const PetscScalar *x_array;
666f650675bSDave May   PetscBool         default_init_guess_value;
667f650675bSDave May 
668f650675bSDave May   PetscFunctionBegin;
669f650675bSDave May   xtmp    = sred->xtmp;
670f650675bSDave May   scatter = sred->scatter;
671f650675bSDave May   yred    = sred->yred;
672f650675bSDave May 
673f650675bSDave May   if (its > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCApplyRichardson_Telescope only supports max_it = 1");
674f650675bSDave May   *reason = (PCRichardsonConvergedReason)0;
675f650675bSDave May 
676f650675bSDave May   if (!zeroguess) {
677f650675bSDave May     ierr = PetscInfo(pc,"PCTelescope: Scattering y for non-zero initial guess\n");CHKERRQ(ierr);
678f650675bSDave May     /* pull in vector y->xtmp */
679f650675bSDave May     ierr = VecScatterBegin(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
680f650675bSDave May     ierr = VecScatterEnd(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
681f650675bSDave May 
682bf00f589SPatrick Sanan     /* copy vector entries into xred */
683f650675bSDave May     ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr);
684f650675bSDave May     if (yred) {
685f650675bSDave May       PetscScalar *LA_yred;
686f650675bSDave May       ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr);
687f650675bSDave May       ierr = VecGetArray(yred,&LA_yred);CHKERRQ(ierr);
688f650675bSDave May       for (i=0; i<ed-st; i++) {
689f650675bSDave May         LA_yred[i] = x_array[i];
690f650675bSDave May       }
691f650675bSDave May       ierr = VecRestoreArray(yred,&LA_yred);CHKERRQ(ierr);
692f650675bSDave May     }
693f650675bSDave May     ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr);
694f650675bSDave May   }
695f650675bSDave May 
69657f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
697f650675bSDave May     ierr = KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value);CHKERRQ(ierr);
698f650675bSDave May     if (!zeroguess) ierr = KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE);CHKERRQ(ierr);
699f650675bSDave May   }
700f650675bSDave May 
701f650675bSDave May   ierr = PCApply_Telescope(pc,x,y);CHKERRQ(ierr);
702f650675bSDave May 
70357f12427SDave May   if (PCTelescope_isActiveRank(sred)) {
704f650675bSDave May     ierr = KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value);CHKERRQ(ierr);
705f650675bSDave May   }
706f650675bSDave May 
707f650675bSDave May   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
708f650675bSDave May   *outits = 1;
709f650675bSDave May   PetscFunctionReturn(0);
710f650675bSDave May }
711f650675bSDave May 
7121e07b27eSBarry Smith static PetscErrorCode PCReset_Telescope(PC pc)
7131e07b27eSBarry Smith {
7141e07b27eSBarry Smith   PC_Telescope   sred = (PC_Telescope)pc->data;
7151e07b27eSBarry Smith   PetscErrorCode ierr;
7161e07b27eSBarry Smith 
7171e07b27eSBarry Smith   ierr = ISDestroy(&sred->isin);CHKERRQ(ierr);
7181e07b27eSBarry Smith   ierr = VecScatterDestroy(&sred->scatter);CHKERRQ(ierr);
719e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->xred);CHKERRQ(ierr);
720e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->yred);CHKERRQ(ierr);
721e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->xtmp);CHKERRQ(ierr);
722e3acf2f7SBarry Smith   ierr = MatDestroy(&sred->Bred);CHKERRQ(ierr);
723e3acf2f7SBarry Smith   ierr = KSPReset(sred->ksp);CHKERRQ(ierr);
7241e07b27eSBarry Smith   if (sred->pctelescope_reset_type) {
7251e07b27eSBarry Smith     ierr = sred->pctelescope_reset_type(pc);CHKERRQ(ierr);
7261e07b27eSBarry Smith   }
7271e07b27eSBarry Smith   PetscFunctionReturn(0);
7281e07b27eSBarry Smith }
7291e07b27eSBarry Smith 
7301e07b27eSBarry Smith static PetscErrorCode PCDestroy_Telescope(PC pc)
7311e07b27eSBarry Smith {
7321e07b27eSBarry Smith   PC_Telescope   sred = (PC_Telescope)pc->data;
7331e07b27eSBarry Smith   PetscErrorCode ierr;
7341e07b27eSBarry Smith 
7351e07b27eSBarry Smith   PetscFunctionBegin;
7361e07b27eSBarry Smith   ierr = PCReset_Telescope(pc);CHKERRQ(ierr);
737e3acf2f7SBarry Smith   ierr = KSPDestroy(&sred->ksp);CHKERRQ(ierr);
7381e07b27eSBarry Smith   ierr = PetscSubcommDestroy(&sred->psubcomm);CHKERRQ(ierr);
739e3acf2f7SBarry Smith   ierr = PetscFree(sred->dm_ctx);CHKERRQ(ierr);
740e3acf2f7SBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
7411e07b27eSBarry Smith   PetscFunctionReturn(0);
7421e07b27eSBarry Smith }
7431e07b27eSBarry Smith 
7444416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Telescope(PetscOptionItems *PetscOptionsObject,PC pc)
7451e07b27eSBarry Smith {
7461e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
7471e07b27eSBarry Smith   PetscErrorCode   ierr;
7481e07b27eSBarry Smith   MPI_Comm         comm;
7491e07b27eSBarry Smith   PetscMPIInt      size;
75048a10b22SPatrick Sanan   PetscBool        flg;
75148a10b22SPatrick Sanan   PetscSubcommType subcommtype;
7521e07b27eSBarry Smith 
7531e07b27eSBarry Smith   PetscFunctionBegin;
7541e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
7551e07b27eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
7561e07b27eSBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Telescope options");CHKERRQ(ierr);
75748a10b22SPatrick Sanan   ierr = PetscOptionsEnum("-pc_telescope_subcomm_type","Subcomm type (interlaced or contiguous)","PCTelescopeSetSubcommType",PetscSubcommTypes,(PetscEnum)sred->subcommtype,(PetscEnum*)&subcommtype,&flg);CHKERRQ(ierr);
75848a10b22SPatrick Sanan   if (flg) {
75948a10b22SPatrick Sanan     ierr = PCTelescopeSetSubcommType(pc,subcommtype);CHKERRQ(ierr);
76048a10b22SPatrick Sanan   }
7611e07b27eSBarry Smith   ierr = PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,0);CHKERRQ(ierr);
7621e07b27eSBarry Smith   if (sred->redfactor > size) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size");
7631e07b27eSBarry Smith   ierr = PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,0);CHKERRQ(ierr);
7647c5279cbSDave May   ierr = PetscOptionsBool("-pc_telescope_ignore_kspcomputeoperators","Ignore method used to compute A","PCTelescopeSetIgnoreKSPComputeOperators",sred->ignore_kspcomputeoperators,&sred->ignore_kspcomputeoperators,0);CHKERRQ(ierr);
7658d9f7141SDave May   ierr = PetscOptionsBool("-pc_telescope_use_coarse_dm","Define sub-communicator from the coarse DM","PCTelescopeSetUseCoarseDM",sred->use_coarse_dm,&sred->use_coarse_dm,0);CHKERRQ(ierr);
7661e07b27eSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7671e07b27eSBarry Smith   PetscFunctionReturn(0);
7681e07b27eSBarry Smith }
7691e07b27eSBarry Smith 
7701e07b27eSBarry Smith /* PC simplementation specific API's */
7711e07b27eSBarry Smith 
7721e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp)
7731e07b27eSBarry Smith {
7741e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
775bd49479cSSatish Balay   PetscFunctionBegin;
7761e07b27eSBarry Smith   if (ksp) *ksp = red->ksp;
777bd49479cSSatish Balay   PetscFunctionReturn(0);
7781e07b27eSBarry Smith }
7791e07b27eSBarry Smith 
78048a10b22SPatrick Sanan static PetscErrorCode PCTelescopeGetSubcommType_Telescope(PC pc,PetscSubcommType *subcommtype)
78148a10b22SPatrick Sanan {
78248a10b22SPatrick Sanan   PC_Telescope red = (PC_Telescope)pc->data;
78348a10b22SPatrick Sanan   PetscFunctionBegin;
78448a10b22SPatrick Sanan   if (subcommtype) *subcommtype = red->subcommtype;
78548a10b22SPatrick Sanan   PetscFunctionReturn(0);
78648a10b22SPatrick Sanan }
78748a10b22SPatrick Sanan 
78848a10b22SPatrick Sanan static PetscErrorCode PCTelescopeSetSubcommType_Telescope(PC pc,PetscSubcommType subcommtype)
78948a10b22SPatrick Sanan {
79048a10b22SPatrick Sanan   PC_Telescope     red = (PC_Telescope)pc->data;
79148a10b22SPatrick Sanan 
79248a10b22SPatrick Sanan   PetscFunctionBegin;
79348a10b22SPatrick Sanan   if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"You cannot change the subcommunicator type for PCTelescope after it has been set up.");
79448a10b22SPatrick Sanan   red->subcommtype = subcommtype;
79548a10b22SPatrick Sanan   PetscFunctionReturn(0);
79648a10b22SPatrick Sanan }
79748a10b22SPatrick Sanan 
7981e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact)
7991e07b27eSBarry Smith {
8001e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
801bd49479cSSatish Balay   PetscFunctionBegin;
8021e07b27eSBarry Smith   if (fact) *fact = red->redfactor;
803bd49479cSSatish Balay   PetscFunctionReturn(0);
8041e07b27eSBarry Smith }
8051e07b27eSBarry Smith 
8061e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact)
8071e07b27eSBarry Smith {
8081e07b27eSBarry Smith   PC_Telescope     red = (PC_Telescope)pc->data;
8091e07b27eSBarry Smith   PetscMPIInt      size;
8101e07b27eSBarry Smith   PetscErrorCode   ierr;
8111e07b27eSBarry Smith 
812bd49479cSSatish Balay   PetscFunctionBegin;
8131e07b27eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
8141e07b27eSBarry Smith   if (fact <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact);
8151e07b27eSBarry Smith   if (fact > size) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact);
8161e07b27eSBarry Smith   red->redfactor = fact;
817bd49479cSSatish Balay   PetscFunctionReturn(0);
8181e07b27eSBarry Smith }
8191e07b27eSBarry Smith 
8201e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v)
8211e07b27eSBarry Smith {
8221e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
823bd49479cSSatish Balay   PetscFunctionBegin;
8241e07b27eSBarry Smith   if (v) *v = red->ignore_dm;
825bd49479cSSatish Balay   PetscFunctionReturn(0);
8261e07b27eSBarry Smith }
82748a10b22SPatrick Sanan 
8281e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v)
8291e07b27eSBarry Smith {
8301e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
831bd49479cSSatish Balay   PetscFunctionBegin;
8321e07b27eSBarry Smith   red->ignore_dm = v;
833bd49479cSSatish Balay   PetscFunctionReturn(0);
8341e07b27eSBarry Smith }
8351e07b27eSBarry Smith 
8368d9f7141SDave May static PetscErrorCode PCTelescopeGetUseCoarseDM_Telescope(PC pc,PetscBool *v)
8378d9f7141SDave May {
8388d9f7141SDave May   PC_Telescope red = (PC_Telescope)pc->data;
8398d9f7141SDave May   PetscFunctionBegin;
8408d9f7141SDave May   if (v) *v = red->use_coarse_dm;
8418d9f7141SDave May   PetscFunctionReturn(0);
8428d9f7141SDave May }
8438d9f7141SDave May 
8448d9f7141SDave May static PetscErrorCode PCTelescopeSetUseCoarseDM_Telescope(PC pc,PetscBool v)
8458d9f7141SDave May {
8468d9f7141SDave May   PC_Telescope red = (PC_Telescope)pc->data;
8478d9f7141SDave May   PetscFunctionBegin;
8488d9f7141SDave May   red->use_coarse_dm = v;
8498d9f7141SDave May   PetscFunctionReturn(0);
8508d9f7141SDave May }
8518d9f7141SDave May 
8520ae7c45bSDave May static PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool *v)
8530ae7c45bSDave May {
8540ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
8550ae7c45bSDave May   PetscFunctionBegin;
8560ae7c45bSDave May   if (v) *v = red->ignore_kspcomputeoperators;
8570ae7c45bSDave May   PetscFunctionReturn(0);
8580ae7c45bSDave May }
85948a10b22SPatrick Sanan 
8600ae7c45bSDave May static PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool v)
8610ae7c45bSDave May {
8620ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
8630ae7c45bSDave May   PetscFunctionBegin;
8640ae7c45bSDave May   red->ignore_kspcomputeoperators = v;
8650ae7c45bSDave May   PetscFunctionReturn(0);
8660ae7c45bSDave May }
8670ae7c45bSDave May 
8681e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm)
8691e07b27eSBarry Smith {
8701e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
871bd49479cSSatish Balay   PetscFunctionBegin;
8721e07b27eSBarry Smith   *dm = private_PCTelescopeGetSubDM(red);
873bd49479cSSatish Balay   PetscFunctionReturn(0);
8741e07b27eSBarry Smith }
8751e07b27eSBarry Smith 
8761e07b27eSBarry Smith /*@
8771e07b27eSBarry Smith  PCTelescopeGetKSP - Gets the KSP created by the telescoping PC.
8781e07b27eSBarry Smith 
8791e07b27eSBarry Smith  Not Collective
8801e07b27eSBarry Smith 
8811e07b27eSBarry Smith  Input Parameter:
8821e07b27eSBarry Smith .  pc - the preconditioner context
8831e07b27eSBarry Smith 
8841e07b27eSBarry Smith  Output Parameter:
8851e07b27eSBarry Smith .  subksp - the KSP defined the smaller set of processes
8861e07b27eSBarry Smith 
8871e07b27eSBarry Smith  Level: advanced
8881e07b27eSBarry Smith 
8891e07b27eSBarry Smith .keywords: PC, telescopting solve
8901e07b27eSBarry Smith @*/
8911e07b27eSBarry Smith PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp)
8921e07b27eSBarry Smith {
893bd49479cSSatish Balay   PetscErrorCode ierr;
894bd49479cSSatish Balay   PetscFunctionBegin;
895163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp));CHKERRQ(ierr);
896bd49479cSSatish Balay   PetscFunctionReturn(0);
8971e07b27eSBarry Smith }
8981e07b27eSBarry Smith 
8991e07b27eSBarry Smith /*@
9001e07b27eSBarry Smith  PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by.
9011e07b27eSBarry Smith 
9021e07b27eSBarry Smith  Not Collective
9031e07b27eSBarry Smith 
9041e07b27eSBarry Smith  Input Parameter:
9051e07b27eSBarry Smith .  pc - the preconditioner context
9061e07b27eSBarry Smith 
9071e07b27eSBarry Smith  Output Parameter:
9081e07b27eSBarry Smith .  fact - the reduction factor
9091e07b27eSBarry Smith 
9101e07b27eSBarry Smith  Level: advanced
9111e07b27eSBarry Smith 
9121e07b27eSBarry Smith .keywords: PC, telescoping solve
9131e07b27eSBarry Smith @*/
9141e07b27eSBarry Smith PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact)
9151e07b27eSBarry Smith {
916bd49479cSSatish Balay   PetscErrorCode ierr;
917bd49479cSSatish Balay   PetscFunctionBegin;
918163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact));CHKERRQ(ierr);
919bd49479cSSatish Balay   PetscFunctionReturn(0);
9201e07b27eSBarry Smith }
9211e07b27eSBarry Smith 
9221e07b27eSBarry Smith /*@
9231e07b27eSBarry Smith  PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by.
9241e07b27eSBarry Smith 
9251e07b27eSBarry Smith  Not Collective
9261e07b27eSBarry Smith 
9271e07b27eSBarry Smith  Input Parameter:
9281e07b27eSBarry Smith .  pc - the preconditioner context
9291e07b27eSBarry Smith 
9301e07b27eSBarry Smith  Output Parameter:
9311e07b27eSBarry Smith .  fact - the reduction factor
9321e07b27eSBarry Smith 
9331e07b27eSBarry Smith  Level: advanced
9341e07b27eSBarry Smith 
9351e07b27eSBarry Smith .keywords: PC, telescoping solve
9361e07b27eSBarry Smith @*/
9371e07b27eSBarry Smith PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact)
9381e07b27eSBarry Smith {
939bd49479cSSatish Balay   PetscErrorCode ierr;
940bd49479cSSatish Balay   PetscFunctionBegin;
941bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact));CHKERRQ(ierr);
942bd49479cSSatish Balay   PetscFunctionReturn(0);
9431e07b27eSBarry Smith }
9441e07b27eSBarry Smith 
9451e07b27eSBarry Smith /*@
9461e07b27eSBarry Smith  PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used.
9471e07b27eSBarry Smith 
9481e07b27eSBarry Smith  Not Collective
9491e07b27eSBarry Smith 
9501e07b27eSBarry Smith  Input Parameter:
9511e07b27eSBarry Smith .  pc - the preconditioner context
9521e07b27eSBarry Smith 
9531e07b27eSBarry Smith  Output Parameter:
9541e07b27eSBarry Smith .  v - the flag
9551e07b27eSBarry Smith 
9561e07b27eSBarry Smith  Level: advanced
9571e07b27eSBarry Smith 
9581e07b27eSBarry Smith .keywords: PC, telescoping solve
9591e07b27eSBarry Smith @*/
9601e07b27eSBarry Smith PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v)
9611e07b27eSBarry Smith {
962bd49479cSSatish Balay   PetscErrorCode ierr;
963bd49479cSSatish Balay   PetscFunctionBegin;
964163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr);
965bd49479cSSatish Balay   PetscFunctionReturn(0);
9661e07b27eSBarry Smith }
9671e07b27eSBarry Smith 
9681e07b27eSBarry Smith /*@
9691e07b27eSBarry Smith  PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC.
9701e07b27eSBarry Smith 
9711e07b27eSBarry Smith  Not Collective
9721e07b27eSBarry Smith 
9731e07b27eSBarry Smith  Input Parameter:
9741e07b27eSBarry Smith .  pc - the preconditioner context
9751e07b27eSBarry Smith 
9761e07b27eSBarry Smith  Output Parameter:
9771e07b27eSBarry Smith .  v - Use PETSC_TRUE to ignore any DM
9781e07b27eSBarry Smith 
9791e07b27eSBarry Smith  Level: advanced
9801e07b27eSBarry Smith 
9811e07b27eSBarry Smith .keywords: PC, telescoping solve
9821e07b27eSBarry Smith @*/
983bfd6bcc6SSatish Balay PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v)
9841e07b27eSBarry Smith {
985bd49479cSSatish Balay   PetscErrorCode ierr;
986bd49479cSSatish Balay   PetscFunctionBegin;
987bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr);
988bd49479cSSatish Balay   PetscFunctionReturn(0);
9891e07b27eSBarry Smith }
9901e07b27eSBarry Smith 
9911e07b27eSBarry Smith /*@
9928d9f7141SDave May  PCTelescopeGetUseCoarseDM - Get the flag indicating if the coarse DM attached to DM associated with the PC will be used.
9938d9f7141SDave May 
9948d9f7141SDave May  Not Collective
9958d9f7141SDave May 
9968d9f7141SDave May  Input Parameter:
9978d9f7141SDave May .  pc - the preconditioner context
9988d9f7141SDave May 
9998d9f7141SDave May  Output Parameter:
10008d9f7141SDave May .  v - the flag
10018d9f7141SDave May 
10028d9f7141SDave May  Level: advanced
10038d9f7141SDave May 
10048d9f7141SDave May .keywords: PC, telescoping solve
10058d9f7141SDave May @*/
10068d9f7141SDave May PetscErrorCode PCTelescopeGetUseCoarseDM(PC pc,PetscBool *v)
10078d9f7141SDave May {
10088d9f7141SDave May   PetscErrorCode ierr;
10098d9f7141SDave May   PetscFunctionBegin;
10108d9f7141SDave May   ierr = PetscUseMethod(pc,"PCTelescopeGetUseCoarseDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr);
10118d9f7141SDave May   PetscFunctionReturn(0);
10128d9f7141SDave May }
10138d9f7141SDave May 
10148d9f7141SDave May /*@
10158d9f7141SDave May  PCTelescopeSetUseCoarseDM - Set a flag to query the DM attached to the PC if it also has a coarse DM
10168d9f7141SDave May 
10178d9f7141SDave May  Not Collective
10188d9f7141SDave May 
10198d9f7141SDave May  Input Parameter:
10208d9f7141SDave May .  pc - the preconditioner context
10218d9f7141SDave May 
10228d9f7141SDave May  Output Parameter:
1023aaa7a805SDave May .  v - Use PETSC_FALSE to ignore any coarse DM
10248d9f7141SDave May 
10258d9f7141SDave May  Notes:
10268d9f7141SDave May  When you have specified to use a coarse DM, the communicator used to create the sub-KSP within PCTelescope
10278d9f7141SDave May  will be that of the coarse DM. Hence the flags -pc_telescope_reduction_factor and
10288d9f7141SDave May  -pc_telescope_subcomm_type will no longer have any meaning.
10298d9f7141SDave May  It is required that the communicator associated with the parent (fine) and the coarse DM are of different sizes.
10308d9f7141SDave May  An error will occur of the size of the communicator associated with the coarse DM
10318d9f7141SDave May  is the same as that of the parent DM.
10328d9f7141SDave May  Furthermore, it is required that the communicator on the coarse DM is a sub-communicator of the parent.
10338d9f7141SDave May  This will be checked at the time the preconditioner is setup and an error will occur if
10348d9f7141SDave May  the coarse DM does not define a sub-communicator of that used by the parent DM.
10358d9f7141SDave May 
10368d9f7141SDave May  The particular Telescope setup invoked when using a coarse DM is agnostic with respect to the type of
10378d9f7141SDave May  the DM used (e.g. it supports DMSHELL, DMPLEX, etc).
10388d9f7141SDave May 
10398d9f7141SDave May  Support is currently only provided for the case when you are using KSPSetComputeOperators()
10408d9f7141SDave May 
10418d9f7141SDave May  The user is required to compose a function with the parent DM to facilitate the transfer of fields (Vec) between the different decompositions defined by the fine and coarse DMs.
10428d9f7141SDave May  In the user code, this is achieved via
10438d9f7141SDave May .vb
10448d9f7141SDave May    {
10458d9f7141SDave May      DM dm_fine;
10468d9f7141SDave May      PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeFieldScatter",your_field_scatter_method);
10478d9f7141SDave May    }
10488d9f7141SDave May .ve
10498d9f7141SDave May  The signature of the user provided field scatter method is
10508d9f7141SDave May .vb
10518d9f7141SDave May    PetscErrorCode your_field_scatter_method(DM dm_fine,Vec x_fine,ScatterMode mode,DM dm_coarse,Vec x_coarse);
10528d9f7141SDave May .ve
10538d9f7141SDave May  The user must provide support for both mode = SCATTER_FORWARD and mode = SCATTER_REVERSE.
10548d9f7141SDave May  SCATTER_FORWARD implies the direction of transfer is from the parent (fine) DM to the coarse DM.
10558d9f7141SDave May 
10568d9f7141SDave May  Optionally, the user may also compose a function with the parent DM to facilitate the transfer
10578d9f7141SDave May  of state variables between the fine and coarse DMs.
10588d9f7141SDave May  In the context of a finite element discretization, an example state variable might be
10598d9f7141SDave May  values associated with quadrature points within each element.
10608d9f7141SDave May  A user provided state scatter method is composed via
10618d9f7141SDave May .vb
10628d9f7141SDave May    {
10638d9f7141SDave May      DM dm_fine;
10648d9f7141SDave May      PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeStateScatter",your_state_scatter_method);
10658d9f7141SDave May    }
10668d9f7141SDave May .ve
10678d9f7141SDave May  The signature of the user provided state scatter method is
10688d9f7141SDave May .vb
10698d9f7141SDave May    PetscErrorCode your_state_scatter_method(DM dm_fine,ScatterMode mode,DM dm_coarse);
10708d9f7141SDave May .ve
10718d9f7141SDave May  SCATTER_FORWARD implies the direction of transfer is from the fine DM to the coarse DM.
10728d9f7141SDave May  The user is only required to support mode = SCATTER_FORWARD.
10738d9f7141SDave May  No assumption is made about the data type of the state variables.
10748d9f7141SDave May  These must be managed by the user and must be accessible from the DM.
10758d9f7141SDave May 
10768d9f7141SDave May  Care must be taken in defining the user context passed to KSPSetComputeOperators() which is to be
10778d9f7141SDave May  associated with the sub-KSP residing within PCTelescope.
10788d9f7141SDave May  In general, PCTelescope assumes that the context on the fine and coarse DM used with
10798d9f7141SDave May  KSPSetComputeOperators() should be "similar" in type or origin.
10808d9f7141SDave May  Specifically the following rules are used to infer what context on the sub-KSP should be.
10818d9f7141SDave May 
10828d9f7141SDave May  First the contexts from the KSP and the fine and coarse DMs are retrieved.
10838d9f7141SDave May  Note that the special case of a DMSHELL context is queried.
10848d9f7141SDave May 
10858d9f7141SDave May .vb
10868d9f7141SDave May    DMKSPGetComputeOperators(dm_fine,&dmfine_kspfunc,&dmfine_kspctx);
10878d9f7141SDave May    DMGetApplicationContext(dm_fine,&dmfine_appctx);
10888d9f7141SDave May    DMShellGetContext(dm_fine,&dmfine_shellctx);
10898d9f7141SDave May 
10908d9f7141SDave May    DMGetApplicationContext(dm_coarse,&dmcoarse_appctx);
10918d9f7141SDave May    DMShellGetContext(dm_coarse,&dmcoarse_shellctx);
10928d9f7141SDave May .ve
10938d9f7141SDave May 
10948d9f7141SDave May  The following rules are then enforced:
10958d9f7141SDave May 
10968d9f7141SDave May  1. If dmfine_kspctx = NULL, then we provide a NULL pointer as the context for the sub-KSP:
10978d9f7141SDave May  KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,NULL);
10988d9f7141SDave May 
10998d9f7141SDave May  2. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_appctx,
11008d9f7141SDave May  check that dmcoarse_appctx is also non-NULL. If this is true, then:
11018d9f7141SDave May  KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,dmcoarse_appctx);
11028d9f7141SDave May 
11038d9f7141SDave May  3. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_shellctx,
11048d9f7141SDave May  check that dmcoarse_shellctx is also non-NULL. If this is true, then:
11058d9f7141SDave May  KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,dmcoarse_shellctx);
11068d9f7141SDave May 
11078d9f7141SDave May  If neither of the above three tests passed, then PCTelescope cannot safely determine what
11088d9f7141SDave May  context should be provided to KSPSetComputeOperators() for use with the sub-KSP.
11098d9f7141SDave May  In this case, an additional mechanism is provided via a composed function which will return
11108d9f7141SDave May  the actual context to be used. To use this feature you must compose the "getter" function
11118d9f7141SDave May  with the coarse DM, e.g.
11128d9f7141SDave May .vb
11138d9f7141SDave May    {
11148d9f7141SDave May      DM dm_coarse;
11158d9f7141SDave May      PetscObjectCompose((PetscObject)dm_coarse,"PCTelescopeGetCoarseDMKSPContext",your_coarse_context_getter);
11168d9f7141SDave May    }
11178d9f7141SDave May .ve
11188d9f7141SDave May  The signature of the user provided method is
11198d9f7141SDave May .vb
11208d9f7141SDave May    PetscErrorCode your_coarse_context_getter(DM dm_coarse,void **your_kspcontext);
11218d9f7141SDave May .ve
11228d9f7141SDave May 
11238d9f7141SDave May  Level: advanced
11248d9f7141SDave May 
11258d9f7141SDave May .keywords: PC, telescoping solve
11268d9f7141SDave May @*/
11278d9f7141SDave May PetscErrorCode PCTelescopeSetUseCoarseDM(PC pc,PetscBool v)
11288d9f7141SDave May {
11298d9f7141SDave May   PetscErrorCode ierr;
11308d9f7141SDave May   PetscFunctionBegin;
11318d9f7141SDave May   ierr = PetscTryMethod(pc,"PCTelescopeSetUseCoarseDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr);
11328d9f7141SDave May   PetscFunctionReturn(0);
11338d9f7141SDave May }
11348d9f7141SDave May 
11358d9f7141SDave May /*@
11360ae7c45bSDave May  PCTelescopeGetIgnoreKSPComputeOperators - Get the flag indicating if KSPComputeOperators will be used.
11370ae7c45bSDave May 
11380ae7c45bSDave May  Not Collective
11390ae7c45bSDave May 
11400ae7c45bSDave May  Input Parameter:
11410ae7c45bSDave May .  pc - the preconditioner context
11420ae7c45bSDave May 
11430ae7c45bSDave May  Output Parameter:
11440ae7c45bSDave May .  v - the flag
11450ae7c45bSDave May 
11460ae7c45bSDave May  Level: advanced
11470ae7c45bSDave May 
11480ae7c45bSDave May .keywords: PC, telescoping solve
11490ae7c45bSDave May @*/
11500ae7c45bSDave May PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC pc,PetscBool *v)
11510ae7c45bSDave May {
11520ae7c45bSDave May   PetscErrorCode ierr;
11530ae7c45bSDave May   PetscFunctionBegin;
1154163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr);
11550ae7c45bSDave May   PetscFunctionReturn(0);
11560ae7c45bSDave May }
11570ae7c45bSDave May 
11580ae7c45bSDave May /*@
11590ae7c45bSDave May  PCTelescopeSetIgnoreKSPComputeOperators - Set a flag to ignore KSPComputeOperators.
11600ae7c45bSDave May 
11610ae7c45bSDave May  Not Collective
11620ae7c45bSDave May 
11630ae7c45bSDave May  Input Parameter:
11640ae7c45bSDave May .  pc - the preconditioner context
11650ae7c45bSDave May 
11660ae7c45bSDave May  Output Parameter:
1167a954d8f4SDave May .  v - Use PETSC_TRUE to ignore the method (if defined) set via KSPSetComputeOperators on pc
11680ae7c45bSDave May 
11690ae7c45bSDave May  Level: advanced
11700ae7c45bSDave May 
11710ae7c45bSDave May .keywords: PC, telescoping solve
11720ae7c45bSDave May @*/
11730ae7c45bSDave May PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC pc,PetscBool v)
11740ae7c45bSDave May {
11750ae7c45bSDave May   PetscErrorCode ierr;
11760ae7c45bSDave May   PetscFunctionBegin;
11770ae7c45bSDave May   ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr);
11780ae7c45bSDave May   PetscFunctionReturn(0);
11790ae7c45bSDave May }
11800ae7c45bSDave May 
11810ae7c45bSDave May /*@
11821e07b27eSBarry Smith  PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP.
11831e07b27eSBarry Smith 
11841e07b27eSBarry Smith  Not Collective
11851e07b27eSBarry Smith 
11861e07b27eSBarry Smith  Input Parameter:
11871e07b27eSBarry Smith .  pc - the preconditioner context
11881e07b27eSBarry Smith 
11891e07b27eSBarry Smith  Output Parameter:
11901e07b27eSBarry Smith .  subdm - The re-partitioned DM
11911e07b27eSBarry Smith 
11921e07b27eSBarry Smith  Level: advanced
11931e07b27eSBarry Smith 
11941e07b27eSBarry Smith .keywords: PC, telescoping solve
11951e07b27eSBarry Smith @*/
11961e07b27eSBarry Smith PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm)
11971e07b27eSBarry Smith {
1198bd49479cSSatish Balay   PetscErrorCode ierr;
1199bd49479cSSatish Balay   PetscFunctionBegin;
1200163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm));CHKERRQ(ierr);
1201bd49479cSSatish Balay   PetscFunctionReturn(0);
12021e07b27eSBarry Smith }
12031e07b27eSBarry Smith 
120448a10b22SPatrick Sanan /*@
120548a10b22SPatrick Sanan  PCTelescopeSetSubcommType - set subcommunicator type (interlaced or contiguous)
120648a10b22SPatrick Sanan 
120748a10b22SPatrick Sanan  Logically Collective
120848a10b22SPatrick Sanan 
120948a10b22SPatrick Sanan  Input Parameter:
12101dae98e4SBarry Smith +  pc - the preconditioner context
12111dae98e4SBarry Smith -  subcommtype - the subcommunicator type (see PetscSubcommType)
121248a10b22SPatrick Sanan 
121348a10b22SPatrick Sanan  Level: advanced
121448a10b22SPatrick Sanan 
121548a10b22SPatrick Sanan .keywords: PC, telescoping solve
121648a10b22SPatrick Sanan 
121748a10b22SPatrick Sanan .seealso: PetscSubcommType, PetscSubcomm, PCTELESCOPE
121848a10b22SPatrick Sanan @*/
121948a10b22SPatrick Sanan PetscErrorCode PCTelescopeSetSubcommType(PC pc, PetscSubcommType subcommtype)
122048a10b22SPatrick Sanan {
122148a10b22SPatrick Sanan   PetscErrorCode ierr;
122248a10b22SPatrick Sanan   PetscFunctionBegin;
122348a10b22SPatrick Sanan   ierr = PetscTryMethod(pc,"PCTelescopeSetSubcommType_C",(PC,PetscSubcommType),(pc,subcommtype));CHKERRQ(ierr);
122448a10b22SPatrick Sanan   PetscFunctionReturn(0);
122548a10b22SPatrick Sanan }
122648a10b22SPatrick Sanan 
122748a10b22SPatrick Sanan /*@
122848a10b22SPatrick Sanan  PCTelescopeGetSubcommType - Get the subcommunicator type (interlaced or contiguous)
122948a10b22SPatrick Sanan 
123048a10b22SPatrick Sanan  Not Collective
123148a10b22SPatrick Sanan 
123248a10b22SPatrick Sanan  Input Parameter:
123348a10b22SPatrick Sanan .  pc - the preconditioner context
123448a10b22SPatrick Sanan 
123548a10b22SPatrick Sanan  Output Parameter:
123648a10b22SPatrick Sanan .  subcommtype - the subcommunicator type (see PetscSubcommType)
123748a10b22SPatrick Sanan 
123848a10b22SPatrick Sanan  Level: advanced
123948a10b22SPatrick Sanan 
124048a10b22SPatrick Sanan .keywords: PC, telescoping solve
124148a10b22SPatrick Sanan 
12421dae98e4SBarry Smith .seealso: PetscSubcomm, PetscSubcommType, PCTELESCOPE
124348a10b22SPatrick Sanan @*/
12441dae98e4SBarry Smith PetscErrorCode PCTelescopeGetSubcommType(PC pc, PetscSubcommType *subcommtype)
124548a10b22SPatrick Sanan {
124648a10b22SPatrick Sanan   PetscErrorCode ierr;
124748a10b22SPatrick Sanan   PetscFunctionBegin;
124848a10b22SPatrick Sanan   ierr = PetscUseMethod(pc,"PCTelescopeGetSubcommType_C",(PC,PetscSubcommType*),(pc,subcommtype));CHKERRQ(ierr);
124948a10b22SPatrick Sanan   PetscFunctionReturn(0);
125048a10b22SPatrick Sanan }
125148a10b22SPatrick Sanan 
12521e07b27eSBarry Smith /* -------------------------------------------------------------------------------------*/
12531e07b27eSBarry Smith /*MC
125400fea0ebSDave May    PCTELESCOPE - Runs a KSP solver on a sub-communicator. MPI ranks not in the sub-communicator are idle during the solve.
12551e07b27eSBarry Smith 
12561e07b27eSBarry Smith    Options Database:
125700fea0ebSDave 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.
125800fea0ebSDave May .  -pc_telescope_ignore_dm  - flag to indicate whether an attached DM should be ignored.
125900fea0ebSDave May .  -pc_telescope_subcomm_type <interlaced,contiguous> - defines the selection of MPI ranks on the sub-communicator. see PetscSubcomm for more information.
126000fea0ebSDave May .  -pc_telescope_ignore_kspcomputeoperators - flag to indicate whether KSPSetComputeOperators should be used on the sub-KSP.
126100fea0ebSDave May -  -pc_telescope_use_coarse_dm - flag to indicate whether the coarse DM should be used to define the sub-communicator.
12621e07b27eSBarry Smith 
12631e07b27eSBarry Smith    Level: advanced
12641e07b27eSBarry Smith 
12651e07b27eSBarry Smith    Notes:
126600fea0ebSDave May    Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation
126700fea0ebSDave May    creates a child sub-communicator (c') containing fewer MPI ranks than the original parent preconditioner (PC).
12686fc41876SBarry Smith    The preconditioner is deemed telescopic as it only calls KSPSolve() on a single
12698439623fSPatrick Sanan    sub-communicator, in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators.
127000fea0ebSDave May    This means there will be MPI ranks which will be idle during the application of this preconditioner.
127100fea0ebSDave May    Additionally, in comparison with PCREDUNDANT, PCTELESCOPE can utilize an attached DM.
12726fc41876SBarry Smith 
127300fea0ebSDave May    The default type of the sub KSP (the KSP defined on c') is PREONLY.
127400fea0ebSDave May 
127500fea0ebSDave May    There are three setup mechanisms for PCTelescope. Features support by each type are described below.
127600fea0ebSDave May    In the following, we will refer to the operators B and B', these are the Bmat provided to the KSP on the
127700fea0ebSDave May    communicators c and c' respectively.
127800fea0ebSDave May 
127900fea0ebSDave May    [1] Default setup
128000fea0ebSDave May    The sub-communicator c' is created via PetscSubcommCreate().
128100fea0ebSDave May    Explicitly defined nullspace and near nullspace vectors will be propogated from B to B'.
128200fea0ebSDave May    Currently there is no support define nullspaces via a user supplied method (e.g. as passed to MatNullSpaceSetFunction()).
128300fea0ebSDave May    No support is provided for KSPSetComputeOperators().
128400fea0ebSDave May    Currently there is no support for the flag -pc_use_amat.
128500fea0ebSDave May 
128600fea0ebSDave May    [2] DM aware setup
128700fea0ebSDave May    If a DM is attached to the PC, it is re-partitioned on the sub-communicator c'.
128800fea0ebSDave May    c' is created via PetscSubcommCreate().
12891e07b27eSBarry Smith    Both the Bmat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned DM.
12901e07b27eSBarry Smith    Currently only support for re-partitioning a DMDA is provided.
129100fea0ebSDave May    Any explicitly defined nullspace or near nullspace vectors attached to the original Bmat operator (B) are extracted, re-partitioned and set on the re-partitioned Bmat operator (B').
129200fea0ebSDave May    Currently there is no support define nullspaces via a user supplied method (e.g. as passed to MatNullSpaceSetFunction()).
129300fea0ebSDave May    Support is provided for KSPSetComputeOperators(). The user provided function and context is propagated to the sub KSP.
129400fea0ebSDave May    This is fragile since the user must ensure that their user context is valid for use on c'.
129500fea0ebSDave May    Currently there is no support for the flag -pc_use_amat.
12961e07b27eSBarry Smith 
129700fea0ebSDave May    [3] Coarse DM setup
129800fea0ebSDave May    If a DM (dmfine) is attached to the PC, dmfine is queried for a "coarse" DM (call this dmcoarse) via DMGetCoarseDM().
129900fea0ebSDave May    PCTELESCOPE will interpret the coarse DM as being defined on a sub-communicator of c.
130000fea0ebSDave May    The communicator associated with dmcoarse will define the c' to be used within PCTELESCOPE.
130100fea0ebSDave May    PCTELESCOPE will check that c' is in fact a sub-communicator of c. If it is not, an error will be reported.
130200fea0ebSDave May    The intention of this setup type is that PCTELESCOPE will use an existing (e.g. user defined) communicator hierarchy, say as would be
130300fea0ebSDave May    available with using multi-grid on unstructured meshes.
130400fea0ebSDave May    This setup will not use the command line options -pc_telescope_reduction_factor or -pc_telescope_subcomm_type.
130500fea0ebSDave May    Any explicitly defined nullspace or near nullspace vectors attached to the original Bmat operator (B) are extracted, scattered into the correct ordering consistent with dmcoarse and set on B'.
130600fea0ebSDave May    Currently there is no support define nullspaces via a user supplied method (e.g. as passed to MatNullSpaceSetFunction()).
130700fea0ebSDave May    There is no general method to permute field orderings, hence only KSPSetComputeOperators() is supported.
130800fea0ebSDave May    The user must use PetscObjectComposeFunction() with dmfine to define the method to scatter fields from dmfine to dmcoarse.
130900fea0ebSDave May    Propogation of the user context for KSPSetComputeOperators() on the sub KSP is attempted by querying the DM contexts associated with dmfine and dmcoarse. Alternatively, the user may use PetscObjectComposeFunction() with dmcoarse to define a method which will return the appropriate user context for KSPSetComputeOperators().
131000fea0ebSDave May    Currently there is no support for the flag -pc_use_amat.
131100fea0ebSDave May    This setup can be invoked by the option -pc_telescope_use_coarse_dm or by calling PCTelescopeSetUseCoarseDM(pc,PETSC_TRUE);
131200fea0ebSDave May    Further information about the user-provided methods required by this setup type are described here PCTelescopeSetUseCoarseDM().
13136fc41876SBarry Smith 
13146fc41876SBarry Smith    Developer Notes:
13156fc41876SBarry Smith    During PCSetup, the B operator is scattered onto c'.
13166fc41876SBarry Smith    Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c').
13178439623fSPatrick Sanan    Then, KSPSolve() is executed on the c' communicator.
13186fc41876SBarry Smith 
13196fc41876SBarry Smith    The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED
1320a04a6428SPatrick Sanan    creation routine by default (this can be changed with -pc_telescope_subcomm_type). We run the sub KSP on only the ranks within the communicator which have a color equal to zero.
13216fc41876SBarry Smith 
1322005d9f20SPatrick Sanan    The telescoping preconditioner is aware of nullspaces and near nullspaces which are attached to the B operator.
13238439623fSPatrick Sanan    In the case where B has a (near) nullspace attached, the (near) nullspace vectors are extracted from B and mapped into
1324005d9f20SPatrick Sanan    a new (near) nullspace, defined on the sub-communicator, which is attached to B' (the B operator which was scattered to c')
13256fc41876SBarry Smith 
132600fea0ebSDave May    The telescoping preconditioner can re-partition an attached DM if it is a DMDA (2D or 3D -
132700fea0ebSDave May    support for 1D DMDAs is not provided). If a DMDA is found, a topolgically equivalent DMDA is created on c'
132800fea0ebSDave May    and this new DM is attached the sub KSP. The design of telescope is such that it should be possible to extend support
13298439623fSPatrick Sanan    for re-partitioning other to DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs.
133000fea0ebSDave May    Alternatively, user-provided re-partitioned DMs can be used via -pc_telescope_use_coarse_dm.
13316fc41876SBarry Smith 
133200fea0ebSDave May    With the default setup mode, B' is defined by fusing rows (in order) associated with MPI ranks common to c and c'.
13336fc41876SBarry Smith 
13348439623fSPatrick Sanan    When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permutation of B
13357dae84e0SHong Zhang    into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatCreateSubMatrices(), (iii) fusing the
13366fc41876SBarry Smith    locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat()
13376fc41876SBarry Smith 
13388439623fSPatrick Sanan    Limitations/improvements include the following.
13398439623fSPatrick Sanan    VecPlaceArray() could be used within PCApply() to improve efficiency and reduce memory usage.
134000fea0ebSDave May    A unified mechanism to query for user contexts as required by KSPSetComputeOperators() and MatNullSpaceSetFunction().
13416fc41876SBarry Smith 
13426fc41876SBarry Smith    The symmetric permutation used when a DMDA is encountered is performed via explicitly assmbleming a permutation matrix P,
13438439623fSPatrick Sanan    and performing P^T.A.P. Possibly it might be more efficient to use MatPermute(). We opted to use P^T.A.P as it appears
13448439623fSPatrick Sanan    VecPermute() does not supported for the use case required here. By computing P, one can permute both the operator and RHS in a
13456fc41876SBarry Smith    consistent manner.
13466fc41876SBarry Smith 
134700fea0ebSDave May    Mapping of vectors (default setup mode) is performed in the following way.
134800fea0ebSDave 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.
13498439623fSPatrick Sanan    Using the interlaced creation routine, the ranks in c with color = 0 will be rank 0 and 2.
135000fea0ebSDave May    We perform the scatter to the sub-communicator in the following way.
1351514db83aSDave May    [1] Given a vector x defined on communicator c
13526fc41876SBarry Smith 
1353514db83aSDave May .vb
1354514db83aSDave May    rank(c)  local values of x
1355514db83aSDave May    ------- ----------------------------------------
1356514db83aSDave May         0   [  0.0,  1.0,  2.0,  3.0,  4.0,  5.0 ]
1357514db83aSDave May         1   [  6.0,  7.0,  8.0,  9.0, 10.0, 11.0 ]
1358514db83aSDave May         2   [ 12.0, 13.0, 14.0, 15.0, 16.0, 17.0 ]
1359514db83aSDave May         3   [ 18.0, 19.0, 20.0, 21.0, 22.0, 23.0 ]
1360514db83aSDave May .ve
13616fc41876SBarry Smith 
1362514db83aSDave May    scatter into xtmp defined also on comm c, so that we have the following values
13636fc41876SBarry Smith 
1364514db83aSDave May .vb
1365514db83aSDave May    rank(c)  local values of xtmp
1366514db83aSDave May    ------- ----------------------------------------------------------------------------
1367514db83aSDave 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 ]
1368514db83aSDave May         1   [ ]
1369514db83aSDave 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 ]
1370514db83aSDave May         3   [ ]
1371514db83aSDave May .ve
13726fc41876SBarry Smith 
13736fc41876SBarry Smith    The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') have no values
13746fc41876SBarry Smith 
13756fc41876SBarry Smith 
1376514db83aSDave 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'.
13776fc41876SBarry Smith    Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0.
13786fc41876SBarry Smith 
1379514db83aSDave May .vb
1380514db83aSDave May    rank(c')  local values of xred
1381514db83aSDave May    -------- ----------------------------------------------------------------------------
1382514db83aSDave 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 ]
1383514db83aSDave 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 ]
1384514db83aSDave May .ve
13851e07b27eSBarry Smith 
13861e07b27eSBarry Smith   Contributed by Dave May
13871e07b27eSBarry Smith 
1388bf00f589SPatrick Sanan   Reference:
1389bf00f589SPatrick Sanan   Dave A. May, Patrick Sanan, Karl Rupp, Matthew G. Knepley, and Barry F. Smith, "Extreme-Scale Multigrid Components within PETSc". 2016. In Proceedings of the Platform for Advanced Scientific Computing Conference (PASC '16). DOI: 10.1145/2929908.2929913
1390bf00f589SPatrick Sanan 
13916fc41876SBarry Smith .seealso:  PCTelescopeGetKSP(), PCTelescopeGetDM(), PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(), PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT
13921e07b27eSBarry Smith M*/
13931e07b27eSBarry Smith PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc)
13941e07b27eSBarry Smith {
13951e07b27eSBarry Smith   PetscErrorCode       ierr;
13961e07b27eSBarry Smith   struct _PC_Telescope *sred;
13971e07b27eSBarry Smith 
13981e07b27eSBarry Smith   PetscFunctionBegin;
13991e07b27eSBarry Smith   ierr = PetscNewLog(pc,&sred);CHKERRQ(ierr);
14002a22aa42SDave May   sred->psubcomm       = NULL;
140148a10b22SPatrick Sanan   sred->subcommtype    = PETSC_SUBCOMM_INTERLACED;
14022a22aa42SDave May   sred->subcomm        = MPI_COMM_NULL;
14031e07b27eSBarry Smith   sred->redfactor      = 1;
14041e07b27eSBarry Smith   sred->ignore_dm      = PETSC_FALSE;
14057c5279cbSDave May   sred->ignore_kspcomputeoperators = PETSC_FALSE;
14068d9f7141SDave May   sred->use_coarse_dm  = PETSC_FALSE;
14071e07b27eSBarry Smith   pc->data             = (void*)sred;
14081e07b27eSBarry Smith 
14091e07b27eSBarry Smith   pc->ops->apply           = PCApply_Telescope;
14101e07b27eSBarry Smith   pc->ops->applytranspose  = NULL;
1411f650675bSDave May   pc->ops->applyrichardson = PCApplyRichardson_Telescope;
14121e07b27eSBarry Smith   pc->ops->setup           = PCSetUp_Telescope;
14131e07b27eSBarry Smith   pc->ops->destroy         = PCDestroy_Telescope;
14141e07b27eSBarry Smith   pc->ops->reset           = PCReset_Telescope;
14151e07b27eSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Telescope;
14161e07b27eSBarry Smith   pc->ops->view            = PCView_Telescope;
14171e07b27eSBarry Smith 
14181e07b27eSBarry Smith   sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
14191e07b27eSBarry Smith   sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
14201e07b27eSBarry Smith   sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
14211e07b27eSBarry Smith   sred->pctelescope_reset_type              = NULL;
14221e07b27eSBarry Smith 
14231e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetKSP_C",PCTelescopeGetKSP_Telescope);CHKERRQ(ierr);
142448a10b22SPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetSubcommType_C",PCTelescopeGetSubcommType_Telescope);CHKERRQ(ierr);
142548a10b22SPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetSubcommType_C",PCTelescopeSetSubcommType_Telescope);CHKERRQ(ierr);
14261e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetReductionFactor_C",PCTelescopeGetReductionFactor_Telescope);CHKERRQ(ierr);
14271e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetReductionFactor_C",PCTelescopeSetReductionFactor_Telescope);CHKERRQ(ierr);
14281e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreDM_C",PCTelescopeGetIgnoreDM_Telescope);CHKERRQ(ierr);
14291e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreDM_C",PCTelescopeSetIgnoreDM_Telescope);CHKERRQ(ierr);
14300ae7c45bSDave May   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",PCTelescopeGetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr);
14310ae7c45bSDave May   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",PCTelescopeSetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr);
14321e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetDM_C",PCTelescopeGetDM_Telescope);CHKERRQ(ierr);
14338d9f7141SDave May   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetUseCoarseDM_C",PCTelescopeGetUseCoarseDM_Telescope);CHKERRQ(ierr);
14348d9f7141SDave May   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetUseCoarseDM_C",PCTelescopeSetUseCoarseDM_Telescope);CHKERRQ(ierr);
14351e07b27eSBarry Smith   PetscFunctionReturn(0);
14361e07b27eSBarry Smith }
1437