xref: /petsc/src/ksp/pc/impls/telescope/telescope.c (revision c5083d92c4ff246cae69ce4f2fa2928232b0a98b)
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"
21bf00f589SPatrick Sanan "  url       = {http://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 /*
31*c5083d92SDave May  default setup mode
321e07b27eSBarry Smith 
33*c5083d92SDave May  [1a] scatter to (FORWARD)
341e07b27eSBarry Smith  x(comm) -> xtmp(comm)
35*c5083d92SDave May  [1b] local copy (to) ranks with color = 0
361e07b27eSBarry Smith  xred(subcomm) <- xtmp
371e07b27eSBarry Smith 
38*c5083d92SDave May  [2] solve on sub KSP to obtain yred(subcomm)
39*c5083d92SDave May 
40*c5083d92SDave May  [3a] local copy (from) ranks with color = 0
411e07b27eSBarry Smith  yred(subcomm) --> xtmp
42*c5083d92SDave May  [2b] scatter from (REVERSE)
431e07b27eSBarry Smith  xtmp(comm) -> y(comm)
441e07b27eSBarry Smith */
451e07b27eSBarry Smith 
462a22aa42SDave May PetscBool PetscSubComm_isActiveRank(PetscSubcomm scomm)
471e07b27eSBarry Smith {
481e07b27eSBarry Smith   if (scomm->color == 0) { return PETSC_TRUE; }
491e07b27eSBarry Smith   else { return PETSC_FALSE; }
501e07b27eSBarry Smith }
511e07b27eSBarry Smith 
522a22aa42SDave May PetscBool isActiveRank(PC_Telescope sred)
532a22aa42SDave May {
542a22aa42SDave May   if (sred->psubcomm) {
552a22aa42SDave May     return(PetscSubComm_isActiveRank(sred->psubcomm));
562a22aa42SDave May   } else {
572a22aa42SDave May     if (sred->subcomm != MPI_COMM_NULL) { return PETSC_TRUE; }
582a22aa42SDave May     else { return PETSC_FALSE; }
592a22aa42SDave May   }
602a22aa42SDave May }
612a22aa42SDave May 
628d9f7141SDave May /*
638d9f7141SDave May   Collective on MPI_Comm[comm_f]
648d9f7141SDave May   Notes
658d9f7141SDave May    * Using comm_f = MPI_COMM_NULL will result in an error
668d9f7141SDave May    * Using comm_c = MPI_COMM_NULL is valid. If all instances of comm_c are NULL the subcomm is not valid.
678d9f7141SDave May    * If any non NULL comm_c communicator cannot map any of its ranks to comm_f, the subcomm is not valid.
688d9f7141SDave May */
698d9f7141SDave May PetscErrorCode PCTelescopeTestValidSubcomm(MPI_Comm comm_f,MPI_Comm comm_c,PetscBool *isvalid)
708d9f7141SDave May {
718d9f7141SDave May   int valid = 1;
728d9f7141SDave May   MPI_Group group_f,group_c;
738d9f7141SDave May   PetscErrorCode ierr;
748d9f7141SDave May   int errorcode;
758d9f7141SDave May   PetscMPIInt count,k,size_f = 0,size_c = 0,size_c_sum = 0;
768d9f7141SDave May   int *ranks_f = NULL,*ranks_c = NULL;
778d9f7141SDave May 
788d9f7141SDave May   if (comm_f == MPI_COMM_NULL) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"comm_f cannot be MPI_COMM_NULL");
798d9f7141SDave May 
808d9f7141SDave May   ierr = MPI_Comm_group(comm_f,&group_f);CHKERRQ(ierr);
818d9f7141SDave May   if (comm_c != MPI_COMM_NULL) {
828d9f7141SDave May     ierr = MPI_Comm_group(comm_c,&group_c);CHKERRQ(ierr);
838d9f7141SDave May   }
848d9f7141SDave May 
858d9f7141SDave May   ierr = MPI_Comm_size(comm_f,&size_f);CHKERRQ(ierr);
868d9f7141SDave May   if (comm_c != MPI_COMM_NULL) {
878d9f7141SDave May     ierr = MPI_Comm_size(comm_c,&size_c);CHKERRQ(ierr);
888d9f7141SDave May   }
898d9f7141SDave May 
908d9f7141SDave May   /* check not all comm_c's are NULL */
918d9f7141SDave May   size_c_sum = size_c;
928d9f7141SDave May   ierr = MPI_Allreduce(MPI_IN_PLACE,&size_c_sum,1,MPI_INT,MPI_SUM,comm_f);CHKERRQ(ierr);
938d9f7141SDave May   if (size_c_sum == 0) {
948d9f7141SDave May     valid = 0;
958d9f7141SDave May   }
968d9f7141SDave May 
978d9f7141SDave May   /* check we can map at least 1 rank in comm_c to comm_f */
988d9f7141SDave May   ierr = PetscMalloc1(size_f,&ranks_f);CHKERRQ(ierr);
998d9f7141SDave May   ierr = PetscMalloc1(size_c,&ranks_c);CHKERRQ(ierr);
1008d9f7141SDave May   for (k=0; k<size_f; k++) {
1018d9f7141SDave May     ranks_f[k] = MPI_UNDEFINED;
1028d9f7141SDave May   }
1038d9f7141SDave May   for (k=0; k<size_c; k++) {
1048d9f7141SDave May     ranks_c[k] = (int)k;
1058d9f7141SDave May   }
1068d9f7141SDave May 
1078d9f7141SDave May   count = 0;
1088d9f7141SDave May   if (comm_c != MPI_COMM_NULL) {
1098d9f7141SDave May     errorcode = MPI_Group_translate_ranks(group_c,size_c,ranks_c,group_f,ranks_f);
1108d9f7141SDave May     for (k=0; k<size_f; k++) {
1118d9f7141SDave May       if (ranks_f[k] == MPI_UNDEFINED) {
1128d9f7141SDave May         count++;
1138d9f7141SDave May       }
1148d9f7141SDave May     }
1158d9f7141SDave May   }
1168d9f7141SDave May   if (count == size_f) {
1178d9f7141SDave May     valid = 0;
1188d9f7141SDave May   }
1198d9f7141SDave May 
1208d9f7141SDave May   ierr = MPI_Allreduce(MPI_IN_PLACE,&valid,1,MPI_INT,MPI_MIN,comm_f);CHKERRQ(ierr);
1218d9f7141SDave May   if (valid == 1) { *isvalid = PETSC_TRUE; }
1228d9f7141SDave May   else { *isvalid = PETSC_FALSE; }
1238d9f7141SDave May 
1248d9f7141SDave May   ierr = PetscFree(ranks_f);CHKERRQ(ierr);
1258d9f7141SDave May   ierr = PetscFree(ranks_c);CHKERRQ(ierr);
1268d9f7141SDave May   ierr = MPI_Group_free(&group_f);CHKERRQ(ierr);
1278d9f7141SDave May   if (comm_c != MPI_COMM_NULL) {
1288d9f7141SDave May     ierr = MPI_Group_free(&group_c);CHKERRQ(ierr);
1298d9f7141SDave May   }
1308d9f7141SDave May   PetscFunctionReturn(0);
1318d9f7141SDave May }
1328d9f7141SDave May 
1331e07b27eSBarry Smith DM private_PCTelescopeGetSubDM(PC_Telescope sred)
1341e07b27eSBarry Smith {
135c6a0d831SBarry Smith   DM subdm = NULL;
1361e07b27eSBarry Smith 
1372a22aa42SDave May   if (!isActiveRank(sred)) { subdm = NULL; }
1381e07b27eSBarry Smith   else {
1391e07b27eSBarry Smith     switch (sred->sr_type) {
1401e07b27eSBarry Smith     case TELESCOPE_DEFAULT: subdm = NULL;
1411e07b27eSBarry Smith       break;
1421e07b27eSBarry Smith     case TELESCOPE_DMDA:    subdm = ((PC_Telescope_DMDACtx*)sred->dm_ctx)->dmrepart;
1431e07b27eSBarry Smith       break;
1441e07b27eSBarry Smith     case TELESCOPE_DMPLEX:  subdm = NULL;
1451e07b27eSBarry Smith       break;
1468d9f7141SDave May     case TELESCOPE_COARSEDM: if (sred->ksp) { KSPGetDM(sred->ksp,&subdm); }
1478d9f7141SDave May       break;
1481e07b27eSBarry Smith     }
1491e07b27eSBarry Smith   }
1501e07b27eSBarry Smith   return(subdm);
1511e07b27eSBarry Smith }
1521e07b27eSBarry Smith 
1531e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_default(PC pc,PC_Telescope sred)
1541e07b27eSBarry Smith {
1551e07b27eSBarry Smith   PetscErrorCode ierr;
1561e07b27eSBarry Smith   PetscInt       m,M,bs,st,ed;
1571e07b27eSBarry Smith   Vec            x,xred,yred,xtmp;
1581e07b27eSBarry Smith   Mat            B;
1591e07b27eSBarry Smith   MPI_Comm       comm,subcomm;
1601e07b27eSBarry Smith   VecScatter     scatter;
1611e07b27eSBarry Smith   IS             isin;
1621e07b27eSBarry Smith 
1631e07b27eSBarry Smith   PetscFunctionBegin;
1641e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: setup (default)\n");CHKERRQ(ierr);
1651e07b27eSBarry Smith   comm = PetscSubcommParent(sred->psubcomm);
1661e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
1671e07b27eSBarry Smith 
1681e07b27eSBarry Smith   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
1691e07b27eSBarry Smith   ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr);
1701e07b27eSBarry Smith   ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr);
1711e07b27eSBarry Smith   ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr);
1721e07b27eSBarry Smith 
1731e07b27eSBarry Smith   xred = NULL;
1743ac26c5eSBarry Smith   m    = 0;
1752a22aa42SDave May   if (isActiveRank(sred)) {
1761e07b27eSBarry Smith     ierr = VecCreate(subcomm,&xred);CHKERRQ(ierr);
1771e07b27eSBarry Smith     ierr = VecSetSizes(xred,PETSC_DECIDE,M);CHKERRQ(ierr);
1781e07b27eSBarry Smith     ierr = VecSetBlockSize(xred,bs);CHKERRQ(ierr);
1791e07b27eSBarry Smith     ierr = VecSetFromOptions(xred);CHKERRQ(ierr);
180ca43db0aSBarry Smith     ierr = VecGetLocalSize(xred,&m);CHKERRQ(ierr);
1811e07b27eSBarry Smith   }
1821e07b27eSBarry Smith 
1831e07b27eSBarry Smith   yred = NULL;
1842a22aa42SDave May   if (isActiveRank(sred)) {
1851e07b27eSBarry Smith     ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr);
1861e07b27eSBarry Smith   }
1871e07b27eSBarry Smith 
1881e07b27eSBarry Smith   ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr);
1891e07b27eSBarry Smith   ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr);
1901e07b27eSBarry Smith   ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr);
1911e07b27eSBarry Smith   ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr);
1921e07b27eSBarry Smith 
1932a22aa42SDave May   if (isActiveRank(sred)) {
1941e07b27eSBarry Smith     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
1951e07b27eSBarry Smith     ierr = ISCreateStride(comm,(ed-st),st,1,&isin);CHKERRQ(ierr);
1961e07b27eSBarry Smith   } else {
1971e07b27eSBarry Smith     ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr);
1983ac26c5eSBarry Smith     ierr = ISCreateStride(comm,0,st,1,&isin);CHKERRQ(ierr);
1991e07b27eSBarry Smith   }
2001e07b27eSBarry Smith   ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr);
2011e07b27eSBarry Smith 
20235928de7SBarry Smith   ierr = VecScatterCreateWithData(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr);
2031e07b27eSBarry Smith 
2041e07b27eSBarry Smith   sred->isin    = isin;
2051e07b27eSBarry Smith   sred->scatter = scatter;
2061e07b27eSBarry Smith   sred->xred    = xred;
2071e07b27eSBarry Smith   sred->yred    = yred;
2081e07b27eSBarry Smith   sred->xtmp    = xtmp;
2091e07b27eSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
2101e07b27eSBarry Smith   PetscFunctionReturn(0);
2111e07b27eSBarry Smith }
2121e07b27eSBarry Smith 
2131e07b27eSBarry Smith PetscErrorCode PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A)
2141e07b27eSBarry Smith {
2151e07b27eSBarry Smith   PetscErrorCode ierr;
2161e07b27eSBarry Smith   MPI_Comm       comm,subcomm;
2171e07b27eSBarry Smith   Mat            Bred,B;
2181e07b27eSBarry Smith   PetscInt       nr,nc;
2191e07b27eSBarry Smith   IS             isrow,iscol;
2201e07b27eSBarry Smith   Mat            Blocal,*_Blocal;
2211e07b27eSBarry Smith 
2221e07b27eSBarry Smith   PetscFunctionBegin;
2231e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (default)\n");CHKERRQ(ierr);
2241e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
2251e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
2261e07b27eSBarry Smith   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
2271e07b27eSBarry Smith   ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr);
2281e07b27eSBarry Smith   isrow = sred->isin;
2291e07b27eSBarry Smith   ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr);
2307dae84e0SHong Zhang   ierr = MatCreateSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr);
2311e07b27eSBarry Smith   Blocal = *_Blocal;
2321e07b27eSBarry Smith   ierr = PetscFree(_Blocal);CHKERRQ(ierr);
2331e07b27eSBarry Smith   Bred = NULL;
2342a22aa42SDave May   if (isActiveRank(sred)) {
2351e07b27eSBarry Smith     PetscInt mm;
2361e07b27eSBarry Smith 
2371e07b27eSBarry Smith     if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; }
2381e07b27eSBarry Smith 
2391e07b27eSBarry Smith     ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr);
2401e07b27eSBarry Smith     ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr);
2411e07b27eSBarry Smith   }
2421e07b27eSBarry Smith   *A = Bred;
2431e07b27eSBarry Smith   ierr = ISDestroy(&iscol);CHKERRQ(ierr);
2441e07b27eSBarry Smith   ierr = MatDestroy(&Blocal);CHKERRQ(ierr);
2451e07b27eSBarry Smith   PetscFunctionReturn(0);
2461e07b27eSBarry Smith }
2471e07b27eSBarry Smith 
248392968a1SPatrick Sanan static PetscErrorCode PCTelescopeSubNullSpaceCreate_Telescope(PC pc,PC_Telescope sred,MatNullSpace nullspace,MatNullSpace *sub_nullspace)
2491e07b27eSBarry Smith {
2501e07b27eSBarry Smith   PetscErrorCode   ierr;
2511e07b27eSBarry Smith   PetscBool        has_const;
2521e07b27eSBarry Smith   const Vec        *vecs;
253c41e779fSDave May   Vec              *sub_vecs = NULL;
254392968a1SPatrick Sanan   PetscInt         i,k,n = 0;
2551e07b27eSBarry Smith   MPI_Comm         subcomm;
2561e07b27eSBarry Smith 
2571e07b27eSBarry Smith   PetscFunctionBegin;
2581e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
2591e07b27eSBarry Smith   ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr);
2601e07b27eSBarry Smith 
2612a22aa42SDave May   if (isActiveRank(sred)) {
262e3acf2f7SBarry Smith     if (n) {
263e3acf2f7SBarry Smith       ierr = VecDuplicateVecs(sred->xred,n,&sub_vecs);CHKERRQ(ierr);
2641e07b27eSBarry Smith     }
2651e07b27eSBarry Smith   }
2661e07b27eSBarry Smith 
2671e07b27eSBarry Smith   /* copy entries */
2681e07b27eSBarry Smith   for (k=0; k<n; k++) {
2691e07b27eSBarry Smith     const PetscScalar *x_array;
2701e07b27eSBarry Smith     PetscScalar       *LA_sub_vec;
27113c30530SDave May     PetscInt          st,ed;
2721e07b27eSBarry Smith 
2731e07b27eSBarry Smith     /* pull in vector x->xtmp */
2741e07b27eSBarry Smith     ierr = VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2751e07b27eSBarry Smith     ierr = VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
27647856c66SBarry Smith     if (sub_vecs) {
277a04a6428SPatrick Sanan       /* copy vector entries into xred */
2781e07b27eSBarry Smith       ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
279ea2b237eSDave May       if (sub_vecs[k]) {
2801e07b27eSBarry Smith         ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr);
2811e07b27eSBarry Smith         ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
2821e07b27eSBarry Smith         for (i=0; i<ed-st; i++) {
2831e07b27eSBarry Smith           LA_sub_vec[i] = x_array[i];
2841e07b27eSBarry Smith         }
2851e07b27eSBarry Smith         ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
2861e07b27eSBarry Smith       }
2871e07b27eSBarry Smith       ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
2881e07b27eSBarry Smith     }
28947856c66SBarry Smith   }
2901e07b27eSBarry Smith 
2912a22aa42SDave May   if (isActiveRank(sred)) {
292d8b9d5b7SPatrick Sanan     /* create new (near) nullspace for redundant object */
293392968a1SPatrick Sanan     ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,sub_nullspace);CHKERRQ(ierr);
294392968a1SPatrick Sanan     ierr = VecDestroyVecs(n,&sub_vecs);CHKERRQ(ierr);
295d8b9d5b7SPatrick Sanan     if (nullspace->remove) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope");
296d8b9d5b7SPatrick 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");
297d8b9d5b7SPatrick Sanan   }
298392968a1SPatrick Sanan   PetscFunctionReturn(0);
299392968a1SPatrick Sanan }
300392968a1SPatrick Sanan 
301392968a1SPatrick Sanan static PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat)
302392968a1SPatrick Sanan {
303392968a1SPatrick Sanan   PetscErrorCode   ierr;
304392968a1SPatrick Sanan   Mat              B;
305392968a1SPatrick Sanan 
306392968a1SPatrick Sanan   PetscFunctionBegin;
307392968a1SPatrick Sanan   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
308392968a1SPatrick Sanan 
309392968a1SPatrick Sanan   /* Propagate the nullspace if it exists */
310392968a1SPatrick Sanan   {
311392968a1SPatrick Sanan     MatNullSpace nullspace,sub_nullspace;
312392968a1SPatrick Sanan     ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr);
313392968a1SPatrick Sanan     if (nullspace) {
314392968a1SPatrick Sanan       ierr = PetscInfo(pc,"PCTelescope: generating nullspace (default)\n");CHKERRQ(ierr);
315392968a1SPatrick Sanan       ierr = PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nullspace,&sub_nullspace);CHKERRQ(ierr);
3162a22aa42SDave May       if (isActiveRank(sred)) {
317392968a1SPatrick Sanan         ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr);
31841ff1ee9SPatrick Sanan         ierr = MatNullSpaceDestroy(&sub_nullspace);CHKERRQ(ierr);
3191e07b27eSBarry Smith       }
320392968a1SPatrick Sanan     }
321392968a1SPatrick Sanan   }
322392968a1SPatrick Sanan 
323392968a1SPatrick Sanan   /* Propagate the near nullspace if it exists */
324392968a1SPatrick Sanan   {
325392968a1SPatrick Sanan     MatNullSpace nearnullspace,sub_nearnullspace;
326392968a1SPatrick Sanan     ierr = MatGetNearNullSpace(B,&nearnullspace);CHKERRQ(ierr);
327392968a1SPatrick Sanan     if (nearnullspace) {
328392968a1SPatrick Sanan       ierr = PetscInfo(pc,"PCTelescope: generating near nullspace (default)\n");CHKERRQ(ierr);
329392968a1SPatrick Sanan       ierr = PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nearnullspace,&sub_nearnullspace);CHKERRQ(ierr);
3302a22aa42SDave May       if (isActiveRank(sred)) {
331392968a1SPatrick Sanan         ierr = MatSetNearNullSpace(sub_mat,sub_nearnullspace);CHKERRQ(ierr);
332392968a1SPatrick Sanan         ierr = MatNullSpaceDestroy(&sub_nearnullspace);CHKERRQ(ierr);
333392968a1SPatrick Sanan       }
334392968a1SPatrick Sanan     }
335392968a1SPatrick Sanan   }
3361e07b27eSBarry Smith   PetscFunctionReturn(0);
3371e07b27eSBarry Smith }
3381e07b27eSBarry Smith 
3391e07b27eSBarry Smith static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer)
3401e07b27eSBarry Smith {
3411e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
3421e07b27eSBarry Smith   PetscErrorCode   ierr;
3431e07b27eSBarry Smith   PetscBool        iascii,isstring;
3441e07b27eSBarry Smith   PetscViewer      subviewer;
3451e07b27eSBarry Smith 
3461e07b27eSBarry Smith   PetscFunctionBegin;
3471e07b27eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
3481e07b27eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
3491e07b27eSBarry Smith   if (iascii) {
3508d9f7141SDave May     {
3511e07b27eSBarry Smith       MPI_Comm    comm,subcomm;
3521e07b27eSBarry Smith       PetscMPIInt comm_size,subcomm_size;
3538d9f7141SDave May       DM          dm = NULL,subdm = NULL;
3541e07b27eSBarry Smith 
3551e07b27eSBarry Smith       ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
3561e07b27eSBarry Smith       subdm = private_PCTelescopeGetSubDM(sred);
3578d9f7141SDave May 
3588d9f7141SDave May       if (sred->psubcomm) {
3591e07b27eSBarry Smith         comm = PetscSubcommParent(sred->psubcomm);
3601e07b27eSBarry Smith         subcomm = PetscSubcommChild(sred->psubcomm);
3611e07b27eSBarry Smith         ierr = MPI_Comm_size(comm,&comm_size);CHKERRQ(ierr);
3621e07b27eSBarry Smith         ierr = MPI_Comm_size(subcomm,&subcomm_size);CHKERRQ(ierr);
3631e07b27eSBarry Smith 
3648d9f7141SDave May         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3658d9f7141SDave May         ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm: parent comm size reduction factor = %D\n",sred->redfactor);CHKERRQ(ierr);
3668d9f7141SDave May         ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm: parent_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size);CHKERRQ(ierr);
36748a10b22SPatrick Sanan         switch (sred->subcommtype) {
36848a10b22SPatrick Sanan           case PETSC_SUBCOMM_INTERLACED :
3698d9f7141SDave May             ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm: type = interlaced\n",sred->subcommtype);CHKERRQ(ierr);
37048a10b22SPatrick Sanan             break;
37148a10b22SPatrick Sanan           case PETSC_SUBCOMM_CONTIGUOUS :
3728d9f7141SDave May             ierr = PetscViewerASCIIPrintf(viewer,"petsc subcomm type = contiguous\n",sred->subcommtype);CHKERRQ(ierr);
37348a10b22SPatrick Sanan             break;
37448a10b22SPatrick Sanan           default :
37548a10b22SPatrick Sanan             SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"General subcomm type not supported by PCTelescope");
37648a10b22SPatrick Sanan         }
3778d9f7141SDave May         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3788d9f7141SDave May       } else {
3798d9f7141SDave May         ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
3808d9f7141SDave May         subcomm = sred->subcomm;
3818d9f7141SDave May         if (!isActiveRank(sred)) {
3828d9f7141SDave May           subcomm = PETSC_COMM_SELF;
3838d9f7141SDave May         }
3848d9f7141SDave May 
3858d9f7141SDave May         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3868d9f7141SDave May         ierr = PetscViewerASCIIPrintf(viewer,"subcomm: using user provided sub-communicator\n");CHKERRQ(ierr);
3878d9f7141SDave May         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3888d9f7141SDave May       }
3898d9f7141SDave May 
3901e07b27eSBarry Smith       ierr = PetscViewerGetSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
3912a22aa42SDave May       if (isActiveRank(sred)) {
3921e07b27eSBarry Smith         ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr);
3931e07b27eSBarry Smith 
3941e07b27eSBarry Smith         if (dm && sred->ignore_dm) {
395efd4aadfSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"ignoring DM\n");CHKERRQ(ierr);
3961e07b27eSBarry Smith         }
3977c5279cbSDave May         if (sred->ignore_kspcomputeoperators) {
398efd4aadfSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"ignoring KSPComputeOperators\n");CHKERRQ(ierr);
3997c5279cbSDave May         }
4001e07b27eSBarry Smith         switch (sred->sr_type) {
4011e07b27eSBarry Smith         case TELESCOPE_DEFAULT:
4028d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"setup type: default\n");CHKERRQ(ierr);
4031e07b27eSBarry Smith           break;
4041e07b27eSBarry Smith         case TELESCOPE_DMDA:
4058d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"setup type: DMDA auto-repartitioning\n");CHKERRQ(ierr);
4068ef9ca65SPatrick Sanan           ierr = DMView_DA_Short(subdm,subviewer);CHKERRQ(ierr);
4071e07b27eSBarry Smith           break;
4081e07b27eSBarry Smith         case TELESCOPE_DMPLEX:
4098d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"setup type: DMPLEX auto-repartitioning\n");CHKERRQ(ierr);
4101e07b27eSBarry Smith           break;
4118d9f7141SDave May         case TELESCOPE_COARSEDM:
4128d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"setup type: coarse DM\n");CHKERRQ(ierr);
4138d9f7141SDave May           break;
4148d9f7141SDave May         }
4158d9f7141SDave May 
4168d9f7141SDave May         if (dm) {
4178d9f7141SDave May           PetscObject obj = (PetscObject)dm;
4188d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"Parent DM object:");CHKERRQ(ierr);
4198d9f7141SDave May           PetscViewerASCIIUseTabs(subviewer,PETSC_FALSE);
4208d9f7141SDave May           if (obj->type_name) { PetscViewerASCIIPrintf(subviewer," type = %s;",obj->type_name); }
4218d9f7141SDave May           if (obj->name) { PetscViewerASCIIPrintf(subviewer," name = %s;",obj->name); }
4228d9f7141SDave May           if (obj->prefix) { PetscViewerASCIIPrintf(subviewer," prefix = %s",obj->prefix); }
4238d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"\n");CHKERRQ(ierr);
4248d9f7141SDave May           PetscViewerASCIIUseTabs(subviewer,PETSC_TRUE);
4258d9f7141SDave May         } else {
4268d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"Parent DM object: NULL\n");CHKERRQ(ierr);
4278d9f7141SDave May         }
4288d9f7141SDave May         if (subdm) {
4298d9f7141SDave May           PetscObject obj = (PetscObject)subdm;
4308d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"Sub DM object:");CHKERRQ(ierr);
4318d9f7141SDave May           PetscViewerASCIIUseTabs(subviewer,PETSC_FALSE);
4328d9f7141SDave May           if (obj->type_name) { PetscViewerASCIIPrintf(subviewer," type = %s;",obj->type_name); }
4338d9f7141SDave May           if (obj->name) { PetscViewerASCIIPrintf(subviewer," name = %s;",obj->name); }
4348d9f7141SDave May           if (obj->prefix) { PetscViewerASCIIPrintf(subviewer," prefix = %s",obj->prefix); }
4358d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"\n");CHKERRQ(ierr);
4368d9f7141SDave May           PetscViewerASCIIUseTabs(subviewer,PETSC_TRUE);
4378d9f7141SDave May         } else {
4388d9f7141SDave May           ierr = PetscViewerASCIIPrintf(subviewer,"Sub DM object: NULL\n");CHKERRQ(ierr);
4391e07b27eSBarry Smith         }
4401e07b27eSBarry Smith 
4411e07b27eSBarry Smith         ierr = KSPView(sred->ksp,subviewer);CHKERRQ(ierr);
4421e07b27eSBarry Smith         ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr);
4431e07b27eSBarry Smith       }
4441e07b27eSBarry Smith       ierr = PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
4451e07b27eSBarry Smith     }
4461e07b27eSBarry Smith   }
4471e07b27eSBarry Smith   PetscFunctionReturn(0);
4481e07b27eSBarry Smith }
4491e07b27eSBarry Smith 
4501e07b27eSBarry Smith static PetscErrorCode PCSetUp_Telescope(PC pc)
4511e07b27eSBarry Smith {
4521e07b27eSBarry Smith   PC_Telescope      sred = (PC_Telescope)pc->data;
4531e07b27eSBarry Smith   PetscErrorCode    ierr;
454bd49479cSSatish Balay   MPI_Comm          comm,subcomm=0;
4551e07b27eSBarry Smith   PCTelescopeType   sr_type;
4561e07b27eSBarry Smith 
4571e07b27eSBarry Smith   PetscFunctionBegin;
4581e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
4591e07b27eSBarry Smith 
4601e07b27eSBarry Smith   /* Determine type of setup/update */
4611e07b27eSBarry Smith   if (!pc->setupcalled) {
4621e07b27eSBarry Smith     PetscBool has_dm,same;
4631e07b27eSBarry Smith     DM        dm;
4641e07b27eSBarry Smith 
4651e07b27eSBarry Smith     sr_type = TELESCOPE_DEFAULT;
4661e07b27eSBarry Smith     has_dm = PETSC_FALSE;
4671e07b27eSBarry Smith     ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
4681e07b27eSBarry Smith     if (dm) { has_dm = PETSC_TRUE; }
4691e07b27eSBarry Smith     if (has_dm) {
4701e07b27eSBarry Smith       /* check for dmda */
4711e07b27eSBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&same);CHKERRQ(ierr);
4721e07b27eSBarry Smith       if (same) {
4731e07b27eSBarry Smith         ierr = PetscInfo(pc,"PCTelescope: found DMDA\n");CHKERRQ(ierr);
4741e07b27eSBarry Smith         sr_type = TELESCOPE_DMDA;
4751e07b27eSBarry Smith       }
4761e07b27eSBarry Smith       /* check for dmplex */
4771e07b27eSBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same);CHKERRQ(ierr);
4781e07b27eSBarry Smith       if (same) {
479994fe344SLisandro Dalcin         ierr = PetscInfo(pc,"PCTelescope: found DMPLEX\n");CHKERRQ(ierr);
4801e07b27eSBarry Smith         sr_type = TELESCOPE_DMPLEX;
4811e07b27eSBarry Smith       }
4828d9f7141SDave May 
4838d9f7141SDave May       if (sred->use_coarse_dm) {
4848d9f7141SDave May         ierr = PetscInfo(pc,"PCTelescope: using coarse DM\n");CHKERRQ(ierr);
4858d9f7141SDave May         sr_type = TELESCOPE_COARSEDM;
4861e07b27eSBarry Smith       }
4871e07b27eSBarry Smith 
4881e07b27eSBarry Smith       if (sred->ignore_dm) {
4898d9f7141SDave May         ierr = PetscInfo(pc,"PCTelescope: ignoring DM\n");CHKERRQ(ierr);
4901e07b27eSBarry Smith         sr_type = TELESCOPE_DEFAULT;
4911e07b27eSBarry Smith       }
4928d9f7141SDave May     }
4931e07b27eSBarry Smith     sred->sr_type = sr_type;
4941e07b27eSBarry Smith   } else {
4951e07b27eSBarry Smith     sr_type = sred->sr_type;
4961e07b27eSBarry Smith   }
4971e07b27eSBarry Smith 
498d8b9d5b7SPatrick Sanan   /* set function pointers for repartition setup, matrix creation/update, matrix (near) nullspace, and reset functionality */
4991e07b27eSBarry Smith   switch (sr_type) {
5001e07b27eSBarry Smith     case TELESCOPE_DEFAULT:
5011e07b27eSBarry Smith       sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
5021e07b27eSBarry Smith       sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
5031e07b27eSBarry Smith       sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
5041e07b27eSBarry Smith       sred->pctelescope_reset_type              = NULL;
5051e07b27eSBarry Smith       break;
5061e07b27eSBarry Smith     case TELESCOPE_DMDA:
5071e07b27eSBarry Smith       pc->ops->apply                            = PCApply_Telescope_dmda;
508f650675bSDave May       pc->ops->applyrichardson                  = PCApplyRichardson_Telescope_dmda;
5091e07b27eSBarry Smith       sred->pctelescope_setup_type              = PCTelescopeSetUp_dmda;
5101e07b27eSBarry Smith       sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_dmda;
5111e07b27eSBarry Smith       sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda;
5121e07b27eSBarry Smith       sred->pctelescope_reset_type              = PCReset_Telescope_dmda;
5131e07b27eSBarry Smith       break;
514d8b9d5b7SPatrick Sanan     case TELESCOPE_DMPLEX: SETERRQ(comm,PETSC_ERR_SUP,"Support for DMPLEX is currently not available");
5151e07b27eSBarry Smith       break;
5168d9f7141SDave May     case TELESCOPE_COARSEDM:
5178d9f7141SDave May       pc->ops->apply                            = PCApply_Telescope_CoarseDM;
5188d9f7141SDave May       pc->ops->applyrichardson                  = PCApplyRichardson_Telescope_CoarseDM;
5198d9f7141SDave May       sred->pctelescope_setup_type              = PCTelescopeSetUp_CoarseDM;
5208d9f7141SDave May       sred->pctelescope_matcreate_type          = NULL;
5218d9f7141SDave May       sred->pctelescope_matnullspacecreate_type = NULL;/*PCTelescopeMatNullSpaceCreate_CoarseDM;*/
5228d9f7141SDave May       sred->pctelescope_reset_type              = PCReset_Telescope_CoarseDM;
5231e07b27eSBarry Smith       break;
5248d9f7141SDave May     default: SETERRQ(comm,PETSC_ERR_SUP,"Support only provided for: repartitioning an operator; repartitioning a DMDA; or using a coarse DM");
5258d9f7141SDave May       break;
5268d9f7141SDave May   }
5278d9f7141SDave May 
5288d9f7141SDave May   /* subcomm definition */
5298d9f7141SDave May   if (!pc->setupcalled) {
5308d9f7141SDave May     if ((sr_type == TELESCOPE_DEFAULT) || (sr_type == TELESCOPE_DMDA)) {
5318d9f7141SDave May       if (!sred->psubcomm) {
5328d9f7141SDave May         ierr = PetscSubcommCreate(comm,&sred->psubcomm);CHKERRQ(ierr);
5338d9f7141SDave May         ierr = PetscSubcommSetNumber(sred->psubcomm,sred->redfactor);CHKERRQ(ierr);
5348d9f7141SDave May         ierr = PetscSubcommSetType(sred->psubcomm,sred->subcommtype);CHKERRQ(ierr);
5358d9f7141SDave May         ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr);
5368d9f7141SDave May         sred->subcomm = PetscSubcommChild(sred->psubcomm);
5378d9f7141SDave May       }
5388d9f7141SDave May     } else { /* query PC for DM, check communicators */
5398d9f7141SDave May       DM          dm,dm_coarse_partition = NULL;
5408d9f7141SDave May       MPI_Comm    comm_fine,comm_coarse_partition = MPI_COMM_NULL;
5418d9f7141SDave May       PetscMPIInt csize_fine=0,csize_coarse_partition=0,cs[2],csg[2],cnt=0;
5428d9f7141SDave May       PetscBool   isvalidsubcomm;
5438d9f7141SDave May 
5448d9f7141SDave May       ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
5458d9f7141SDave May       comm_fine = PetscObjectComm((PetscObject)dm);
5468d9f7141SDave May       ierr = DMGetCoarseDM(dm,&dm_coarse_partition);CHKERRQ(ierr);
5478d9f7141SDave May       if (dm_coarse_partition) { cnt = 1; }
5488d9f7141SDave May       ierr = MPI_Allreduce(MPI_IN_PLACE,&cnt,1,MPI_INT,MPI_SUM,comm_fine);CHKERRQ(ierr);
5498d9f7141SDave May       if (cnt == 0) SETERRQ(comm_fine,PETSC_ERR_SUP,"Zero instances of a coarse DM were found");
5508d9f7141SDave May 
5518d9f7141SDave May       ierr = MPI_Comm_size(comm_fine,&csize_fine);CHKERRQ(ierr);
5528d9f7141SDave May       if (dm_coarse_partition) {
5538d9f7141SDave May         comm_coarse_partition = PetscObjectComm((PetscObject)dm_coarse_partition);
5548d9f7141SDave May         ierr = MPI_Comm_size(comm_coarse_partition,&csize_coarse_partition);CHKERRQ(ierr);
5558d9f7141SDave May       }
5568d9f7141SDave May 
5578d9f7141SDave May       cs[0] = csize_fine;
5588d9f7141SDave May       cs[1] = csize_coarse_partition;
5598d9f7141SDave May       ierr = MPI_Allreduce(cs,csg,2,MPI_INT,MPI_MAX,comm_fine);CHKERRQ(ierr);
5608d9f7141SDave 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");
5618d9f7141SDave May 
5628d9f7141SDave May       ierr = PCTelescopeTestValidSubcomm(comm_fine,comm_coarse_partition,&isvalidsubcomm);CHKERRQ(ierr);
5638d9f7141SDave May       if (!isvalidsubcomm) SETERRQ(comm_fine,PETSC_ERR_SUP,"Coarse DM communicator is not a sub-communicator of parentDM->comm");
5648d9f7141SDave May       sred->subcomm = comm_coarse_partition;
5658d9f7141SDave May     }
5668d9f7141SDave May   }
5678d9f7141SDave May   subcomm = sred->subcomm;
5688d9f7141SDave May 
5698d9f7141SDave May   /* internal KSP */
5708d9f7141SDave May   if (!pc->setupcalled) {
5718d9f7141SDave May     const char *prefix;
5728d9f7141SDave May 
5738d9f7141SDave May     if (isActiveRank(sred)) {
5748d9f7141SDave May       ierr = KSPCreate(subcomm,&sred->ksp);CHKERRQ(ierr);
5758d9f7141SDave May       ierr = KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure);CHKERRQ(ierr);
5768d9f7141SDave May       ierr = PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
5778d9f7141SDave May       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp);CHKERRQ(ierr);
5788d9f7141SDave May       ierr = KSPSetType(sred->ksp,KSPPREONLY);CHKERRQ(ierr);
5798d9f7141SDave May       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
5808d9f7141SDave May       ierr = KSPSetOptionsPrefix(sred->ksp,prefix);CHKERRQ(ierr);
5818d9f7141SDave May       ierr = KSPAppendOptionsPrefix(sred->ksp,"telescope_");CHKERRQ(ierr);
5828d9f7141SDave May     }
5831e07b27eSBarry Smith   }
5841e07b27eSBarry Smith 
5851e07b27eSBarry Smith   /* setup */
5861e07b27eSBarry Smith   if (sred->pctelescope_setup_type) {
5871e07b27eSBarry Smith     ierr = sred->pctelescope_setup_type(pc,sred);CHKERRQ(ierr);
5881e07b27eSBarry Smith   }
5891e07b27eSBarry Smith   /* update */
5901e07b27eSBarry Smith   if (!pc->setupcalled) {
5911e07b27eSBarry Smith     if (sred->pctelescope_matcreate_type) {
5921e07b27eSBarry Smith       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred);CHKERRQ(ierr);
5931e07b27eSBarry Smith     }
5941e07b27eSBarry Smith     if (sred->pctelescope_matnullspacecreate_type) {
595392968a1SPatrick Sanan       ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr);
5961e07b27eSBarry Smith     }
5971e07b27eSBarry Smith   } else {
5981e07b27eSBarry Smith     if (sred->pctelescope_matcreate_type) {
5991e07b27eSBarry Smith       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred);CHKERRQ(ierr);
6001e07b27eSBarry Smith     }
6011e07b27eSBarry Smith   }
6021e07b27eSBarry Smith 
6031e07b27eSBarry Smith   /* common - no construction */
6042a22aa42SDave May   if (isActiveRank(sred)) {
6051e07b27eSBarry Smith     ierr = KSPSetOperators(sred->ksp,sred->Bred,sred->Bred);CHKERRQ(ierr);
6061e07b27eSBarry Smith     if (pc->setfromoptionscalled && !pc->setupcalled){
6071e07b27eSBarry Smith       ierr = KSPSetFromOptions(sred->ksp);CHKERRQ(ierr);
6081e07b27eSBarry Smith     }
6091e07b27eSBarry Smith   }
6108d9f7141SDave May 
6118d9f7141SDave May #if 0
6128d9f7141SDave May   /* we perform this last as Bred is not available with KSPSetComputeOperators() until KSPSetUp has been called */
6138d9f7141SDave May   if (!pc->setupcalled) {
6148d9f7141SDave May     if (isActiveRank(sred)) {
6158d9f7141SDave May       ierr = KSPSetUp(sred->ksp);CHKERRQ(ierr);
6168d9f7141SDave May     }
6178d9f7141SDave May     if (sred->pctelescope_matnullspacecreate_type) {
6188d9f7141SDave May       ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr);
6198d9f7141SDave May     }
6208d9f7141SDave May   }
6218d9f7141SDave May #endif
6228d9f7141SDave May 
6231e07b27eSBarry Smith   PetscFunctionReturn(0);
6241e07b27eSBarry Smith }
6251e07b27eSBarry Smith 
6261e07b27eSBarry Smith static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y)
6271e07b27eSBarry Smith {
6281e07b27eSBarry Smith   PC_Telescope      sred = (PC_Telescope)pc->data;
6291e07b27eSBarry Smith   PetscErrorCode    ierr;
6301e07b27eSBarry Smith   Vec               xtmp,xred,yred;
63113c30530SDave May   PetscInt          i,st,ed;
6321e07b27eSBarry Smith   VecScatter        scatter;
6331e07b27eSBarry Smith   PetscScalar       *array;
6341e07b27eSBarry Smith   const PetscScalar *x_array;
6351e07b27eSBarry Smith 
6361e07b27eSBarry Smith   PetscFunctionBegin;
637bf00f589SPatrick Sanan   ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr);
638bf00f589SPatrick Sanan 
6391e07b27eSBarry Smith   xtmp    = sred->xtmp;
6401e07b27eSBarry Smith   scatter = sred->scatter;
6411e07b27eSBarry Smith   xred    = sred->xred;
6421e07b27eSBarry Smith   yred    = sred->yred;
6431e07b27eSBarry Smith 
6441e07b27eSBarry Smith   /* pull in vector x->xtmp */
6451e07b27eSBarry Smith   ierr = VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6461e07b27eSBarry Smith   ierr = VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6471e07b27eSBarry Smith 
648bf00f589SPatrick Sanan   /* copy vector entries into xred */
6491e07b27eSBarry Smith   ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr);
6501e07b27eSBarry Smith   if (xred) {
6511e07b27eSBarry Smith     PetscScalar *LA_xred;
6521e07b27eSBarry Smith     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
6531e07b27eSBarry Smith     ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr);
6541e07b27eSBarry Smith     for (i=0; i<ed-st; i++) {
6551e07b27eSBarry Smith       LA_xred[i] = x_array[i];
6561e07b27eSBarry Smith     }
6571e07b27eSBarry Smith     ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr);
6581e07b27eSBarry Smith   }
6591e07b27eSBarry Smith   ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr);
6601e07b27eSBarry Smith   /* solve */
6612a22aa42SDave May   if (isActiveRank(sred)) {
6621e07b27eSBarry Smith     ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr);
663c0decd05SBarry Smith     ierr = KSPCheckSolve(sred->ksp,pc,yred);CHKERRQ(ierr);
6641e07b27eSBarry Smith   }
6651e07b27eSBarry Smith   /* return vector */
6661e07b27eSBarry Smith   ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr);
6671e07b27eSBarry Smith   if (yred) {
6681e07b27eSBarry Smith     const PetscScalar *LA_yred;
6691e07b27eSBarry Smith     ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr);
6701e07b27eSBarry Smith     ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr);
6711e07b27eSBarry Smith     for (i=0; i<ed-st; i++) {
6721e07b27eSBarry Smith       array[i] = LA_yred[i];
6731e07b27eSBarry Smith     }
6741e07b27eSBarry Smith     ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr);
6751e07b27eSBarry Smith   }
6761e07b27eSBarry Smith   ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr);
6771e07b27eSBarry Smith   ierr = VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6781e07b27eSBarry Smith   ierr = VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
6791e07b27eSBarry Smith   PetscFunctionReturn(0);
6801e07b27eSBarry Smith }
6811e07b27eSBarry Smith 
682f650675bSDave 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)
683f650675bSDave May {
684f650675bSDave May   PC_Telescope      sred = (PC_Telescope)pc->data;
685f650675bSDave May   PetscErrorCode    ierr;
686a1d91a28SDave May   Vec               xtmp,yred;
687f650675bSDave May   PetscInt          i,st,ed;
688f650675bSDave May   VecScatter        scatter;
689f650675bSDave May   const PetscScalar *x_array;
690f650675bSDave May   PetscBool         default_init_guess_value;
691f650675bSDave May 
692f650675bSDave May   PetscFunctionBegin;
693f650675bSDave May   xtmp    = sred->xtmp;
694f650675bSDave May   scatter = sred->scatter;
695f650675bSDave May   yred    = sred->yred;
696f650675bSDave May 
697f650675bSDave May   if (its > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCApplyRichardson_Telescope only supports max_it = 1");
698f650675bSDave May   *reason = (PCRichardsonConvergedReason)0;
699f650675bSDave May 
700f650675bSDave May   if (!zeroguess) {
701f650675bSDave May     ierr = PetscInfo(pc,"PCTelescope: Scattering y for non-zero initial guess\n");CHKERRQ(ierr);
702f650675bSDave May     /* pull in vector y->xtmp */
703f650675bSDave May     ierr = VecScatterBegin(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
704f650675bSDave May     ierr = VecScatterEnd(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
705f650675bSDave May 
706bf00f589SPatrick Sanan     /* copy vector entries into xred */
707f650675bSDave May     ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr);
708f650675bSDave May     if (yred) {
709f650675bSDave May       PetscScalar *LA_yred;
710f650675bSDave May       ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr);
711f650675bSDave May       ierr = VecGetArray(yred,&LA_yred);CHKERRQ(ierr);
712f650675bSDave May       for (i=0; i<ed-st; i++) {
713f650675bSDave May         LA_yred[i] = x_array[i];
714f650675bSDave May       }
715f650675bSDave May       ierr = VecRestoreArray(yred,&LA_yred);CHKERRQ(ierr);
716f650675bSDave May     }
717f650675bSDave May     ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr);
718f650675bSDave May   }
719f650675bSDave May 
7202a22aa42SDave May   if (isActiveRank(sred)) {
721f650675bSDave May     ierr = KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value);CHKERRQ(ierr);
722f650675bSDave May     if (!zeroguess) ierr = KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE);CHKERRQ(ierr);
723f650675bSDave May   }
724f650675bSDave May 
725f650675bSDave May   ierr = PCApply_Telescope(pc,x,y);CHKERRQ(ierr);
726f650675bSDave May 
7272a22aa42SDave May   if (isActiveRank(sred)) {
728f650675bSDave May     ierr = KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value);CHKERRQ(ierr);
729f650675bSDave May   }
730f650675bSDave May 
731f650675bSDave May   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
732f650675bSDave May   *outits = 1;
733f650675bSDave May   PetscFunctionReturn(0);
734f650675bSDave May }
735f650675bSDave May 
7361e07b27eSBarry Smith static PetscErrorCode PCReset_Telescope(PC pc)
7371e07b27eSBarry Smith {
7381e07b27eSBarry Smith   PC_Telescope   sred = (PC_Telescope)pc->data;
7391e07b27eSBarry Smith   PetscErrorCode ierr;
7401e07b27eSBarry Smith 
7411e07b27eSBarry Smith   ierr = ISDestroy(&sred->isin);CHKERRQ(ierr);
7421e07b27eSBarry Smith   ierr = VecScatterDestroy(&sred->scatter);CHKERRQ(ierr);
743e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->xred);CHKERRQ(ierr);
744e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->yred);CHKERRQ(ierr);
745e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->xtmp);CHKERRQ(ierr);
746e3acf2f7SBarry Smith   ierr = MatDestroy(&sred->Bred);CHKERRQ(ierr);
747e3acf2f7SBarry Smith   ierr = KSPReset(sred->ksp);CHKERRQ(ierr);
7481e07b27eSBarry Smith   if (sred->pctelescope_reset_type) {
7491e07b27eSBarry Smith     ierr = sred->pctelescope_reset_type(pc);CHKERRQ(ierr);
7501e07b27eSBarry Smith   }
7511e07b27eSBarry Smith   PetscFunctionReturn(0);
7521e07b27eSBarry Smith }
7531e07b27eSBarry Smith 
7541e07b27eSBarry Smith static PetscErrorCode PCDestroy_Telescope(PC pc)
7551e07b27eSBarry Smith {
7561e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
7571e07b27eSBarry Smith   PetscErrorCode   ierr;
7581e07b27eSBarry Smith 
7591e07b27eSBarry Smith   PetscFunctionBegin;
7601e07b27eSBarry Smith   ierr = PCReset_Telescope(pc);CHKERRQ(ierr);
761e3acf2f7SBarry Smith   ierr = KSPDestroy(&sred->ksp);CHKERRQ(ierr);
7621e07b27eSBarry Smith   ierr = PetscSubcommDestroy(&sred->psubcomm);CHKERRQ(ierr);
763e3acf2f7SBarry Smith   ierr = PetscFree(sred->dm_ctx);CHKERRQ(ierr);
764e3acf2f7SBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
7651e07b27eSBarry Smith   PetscFunctionReturn(0);
7661e07b27eSBarry Smith }
7671e07b27eSBarry Smith 
7684416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Telescope(PetscOptionItems *PetscOptionsObject,PC pc)
7691e07b27eSBarry Smith {
7701e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
7711e07b27eSBarry Smith   PetscErrorCode   ierr;
7721e07b27eSBarry Smith   MPI_Comm         comm;
7731e07b27eSBarry Smith   PetscMPIInt      size;
77448a10b22SPatrick Sanan   PetscBool        flg;
77548a10b22SPatrick Sanan   PetscSubcommType subcommtype;
7761e07b27eSBarry Smith 
7771e07b27eSBarry Smith   PetscFunctionBegin;
7781e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
7791e07b27eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
7801e07b27eSBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Telescope options");CHKERRQ(ierr);
78148a10b22SPatrick Sanan   ierr = PetscOptionsEnum("-pc_telescope_subcomm_type","Subcomm type (interlaced or contiguous)","PCTelescopeSetSubcommType",PetscSubcommTypes,(PetscEnum)sred->subcommtype,(PetscEnum*)&subcommtype,&flg);CHKERRQ(ierr);
78248a10b22SPatrick Sanan   if (flg) {
78348a10b22SPatrick Sanan     ierr = PCTelescopeSetSubcommType(pc,subcommtype);CHKERRQ(ierr);
78448a10b22SPatrick Sanan   }
7851e07b27eSBarry Smith   ierr = PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,0);CHKERRQ(ierr);
7861e07b27eSBarry Smith   if (sred->redfactor > size) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size");
7871e07b27eSBarry Smith   ierr = PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,0);CHKERRQ(ierr);
7887c5279cbSDave May   ierr = PetscOptionsBool("-pc_telescope_ignore_kspcomputeoperators","Ignore method used to compute A","PCTelescopeSetIgnoreKSPComputeOperators",sred->ignore_kspcomputeoperators,&sred->ignore_kspcomputeoperators,0);CHKERRQ(ierr);
7898d9f7141SDave 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);
7901e07b27eSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
7911e07b27eSBarry Smith   PetscFunctionReturn(0);
7921e07b27eSBarry Smith }
7931e07b27eSBarry Smith 
7941e07b27eSBarry Smith /* PC simplementation specific API's */
7951e07b27eSBarry Smith 
7961e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp)
7971e07b27eSBarry Smith {
7981e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
799bd49479cSSatish Balay   PetscFunctionBegin;
8001e07b27eSBarry Smith   if (ksp) *ksp = red->ksp;
801bd49479cSSatish Balay   PetscFunctionReturn(0);
8021e07b27eSBarry Smith }
8031e07b27eSBarry Smith 
80448a10b22SPatrick Sanan static PetscErrorCode PCTelescopeGetSubcommType_Telescope(PC pc,PetscSubcommType *subcommtype)
80548a10b22SPatrick Sanan {
80648a10b22SPatrick Sanan   PC_Telescope red = (PC_Telescope)pc->data;
80748a10b22SPatrick Sanan   PetscFunctionBegin;
80848a10b22SPatrick Sanan   if (subcommtype) *subcommtype = red->subcommtype;
80948a10b22SPatrick Sanan   PetscFunctionReturn(0);
81048a10b22SPatrick Sanan }
81148a10b22SPatrick Sanan 
81248a10b22SPatrick Sanan static PetscErrorCode PCTelescopeSetSubcommType_Telescope(PC pc,PetscSubcommType subcommtype)
81348a10b22SPatrick Sanan {
81448a10b22SPatrick Sanan   PC_Telescope     red = (PC_Telescope)pc->data;
81548a10b22SPatrick Sanan 
81648a10b22SPatrick Sanan   PetscFunctionBegin;
81748a10b22SPatrick 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.");
81848a10b22SPatrick Sanan   red->subcommtype = subcommtype;
81948a10b22SPatrick Sanan   PetscFunctionReturn(0);
82048a10b22SPatrick Sanan }
82148a10b22SPatrick Sanan 
8221e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact)
8231e07b27eSBarry Smith {
8241e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
825bd49479cSSatish Balay   PetscFunctionBegin;
8261e07b27eSBarry Smith   if (fact) *fact = red->redfactor;
827bd49479cSSatish Balay   PetscFunctionReturn(0);
8281e07b27eSBarry Smith }
8291e07b27eSBarry Smith 
8301e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact)
8311e07b27eSBarry Smith {
8321e07b27eSBarry Smith   PC_Telescope     red = (PC_Telescope)pc->data;
8331e07b27eSBarry Smith   PetscMPIInt      size;
8341e07b27eSBarry Smith   PetscErrorCode   ierr;
8351e07b27eSBarry Smith 
836bd49479cSSatish Balay   PetscFunctionBegin;
8371e07b27eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
8381e07b27eSBarry Smith   if (fact <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact);
8391e07b27eSBarry Smith   if (fact > size) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact);
8401e07b27eSBarry Smith   red->redfactor = fact;
841bd49479cSSatish Balay   PetscFunctionReturn(0);
8421e07b27eSBarry Smith }
8431e07b27eSBarry Smith 
8441e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v)
8451e07b27eSBarry Smith {
8461e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
847bd49479cSSatish Balay   PetscFunctionBegin;
8481e07b27eSBarry Smith   if (v) *v = red->ignore_dm;
849bd49479cSSatish Balay   PetscFunctionReturn(0);
8501e07b27eSBarry Smith }
85148a10b22SPatrick Sanan 
8521e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v)
8531e07b27eSBarry Smith {
8541e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
855bd49479cSSatish Balay   PetscFunctionBegin;
8561e07b27eSBarry Smith   red->ignore_dm = v;
857bd49479cSSatish Balay   PetscFunctionReturn(0);
8581e07b27eSBarry Smith }
8591e07b27eSBarry Smith 
8608d9f7141SDave May static PetscErrorCode PCTelescopeGetUseCoarseDM_Telescope(PC pc,PetscBool *v)
8618d9f7141SDave May {
8628d9f7141SDave May   PC_Telescope red = (PC_Telescope)pc->data;
8638d9f7141SDave May   PetscFunctionBegin;
8648d9f7141SDave May   if (v) *v = red->use_coarse_dm;
8658d9f7141SDave May   PetscFunctionReturn(0);
8668d9f7141SDave May }
8678d9f7141SDave May 
8688d9f7141SDave May static PetscErrorCode PCTelescopeSetUseCoarseDM_Telescope(PC pc,PetscBool v)
8698d9f7141SDave May {
8708d9f7141SDave May   PC_Telescope red = (PC_Telescope)pc->data;
8718d9f7141SDave May   PetscFunctionBegin;
8728d9f7141SDave May   red->use_coarse_dm = v;
8738d9f7141SDave May   PetscFunctionReturn(0);
8748d9f7141SDave May }
8758d9f7141SDave May 
8760ae7c45bSDave May static PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool *v)
8770ae7c45bSDave May {
8780ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
8790ae7c45bSDave May   PetscFunctionBegin;
8800ae7c45bSDave May   if (v) *v = red->ignore_kspcomputeoperators;
8810ae7c45bSDave May   PetscFunctionReturn(0);
8820ae7c45bSDave May }
88348a10b22SPatrick Sanan 
8840ae7c45bSDave May static PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool v)
8850ae7c45bSDave May {
8860ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
8870ae7c45bSDave May   PetscFunctionBegin;
8880ae7c45bSDave May   red->ignore_kspcomputeoperators = v;
8890ae7c45bSDave May   PetscFunctionReturn(0);
8900ae7c45bSDave May }
8910ae7c45bSDave May 
8921e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm)
8931e07b27eSBarry Smith {
8941e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
895bd49479cSSatish Balay   PetscFunctionBegin;
8961e07b27eSBarry Smith   *dm = private_PCTelescopeGetSubDM(red);
897bd49479cSSatish Balay   PetscFunctionReturn(0);
8981e07b27eSBarry Smith }
8991e07b27eSBarry Smith 
9001e07b27eSBarry Smith /*@
9011e07b27eSBarry Smith  PCTelescopeGetKSP - Gets the KSP created by the telescoping PC.
9021e07b27eSBarry Smith 
9031e07b27eSBarry Smith  Not Collective
9041e07b27eSBarry Smith 
9051e07b27eSBarry Smith  Input Parameter:
9061e07b27eSBarry Smith .  pc - the preconditioner context
9071e07b27eSBarry Smith 
9081e07b27eSBarry Smith  Output Parameter:
9091e07b27eSBarry Smith .  subksp - the KSP defined the smaller set of processes
9101e07b27eSBarry Smith 
9111e07b27eSBarry Smith  Level: advanced
9121e07b27eSBarry Smith 
9131e07b27eSBarry Smith .keywords: PC, telescopting solve
9141e07b27eSBarry Smith @*/
9151e07b27eSBarry Smith PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp)
9161e07b27eSBarry Smith {
917bd49479cSSatish Balay   PetscErrorCode ierr;
918bd49479cSSatish Balay   PetscFunctionBegin;
919163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp));CHKERRQ(ierr);
920bd49479cSSatish Balay   PetscFunctionReturn(0);
9211e07b27eSBarry Smith }
9221e07b27eSBarry Smith 
9231e07b27eSBarry Smith /*@
9241e07b27eSBarry Smith  PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by.
9251e07b27eSBarry Smith 
9261e07b27eSBarry Smith  Not Collective
9271e07b27eSBarry Smith 
9281e07b27eSBarry Smith  Input Parameter:
9291e07b27eSBarry Smith .  pc - the preconditioner context
9301e07b27eSBarry Smith 
9311e07b27eSBarry Smith  Output Parameter:
9321e07b27eSBarry Smith .  fact - the reduction factor
9331e07b27eSBarry Smith 
9341e07b27eSBarry Smith  Level: advanced
9351e07b27eSBarry Smith 
9361e07b27eSBarry Smith .keywords: PC, telescoping solve
9371e07b27eSBarry Smith @*/
9381e07b27eSBarry Smith PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact)
9391e07b27eSBarry Smith {
940bd49479cSSatish Balay   PetscErrorCode ierr;
941bd49479cSSatish Balay   PetscFunctionBegin;
942163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact));CHKERRQ(ierr);
943bd49479cSSatish Balay   PetscFunctionReturn(0);
9441e07b27eSBarry Smith }
9451e07b27eSBarry Smith 
9461e07b27eSBarry Smith /*@
9471e07b27eSBarry Smith  PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by.
9481e07b27eSBarry Smith 
9491e07b27eSBarry Smith  Not Collective
9501e07b27eSBarry Smith 
9511e07b27eSBarry Smith  Input Parameter:
9521e07b27eSBarry Smith .  pc - the preconditioner context
9531e07b27eSBarry Smith 
9541e07b27eSBarry Smith  Output Parameter:
9551e07b27eSBarry Smith .  fact - the reduction factor
9561e07b27eSBarry Smith 
9571e07b27eSBarry Smith  Level: advanced
9581e07b27eSBarry Smith 
9591e07b27eSBarry Smith .keywords: PC, telescoping solve
9601e07b27eSBarry Smith @*/
9611e07b27eSBarry Smith PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact)
9621e07b27eSBarry Smith {
963bd49479cSSatish Balay   PetscErrorCode ierr;
964bd49479cSSatish Balay   PetscFunctionBegin;
965bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact));CHKERRQ(ierr);
966bd49479cSSatish Balay   PetscFunctionReturn(0);
9671e07b27eSBarry Smith }
9681e07b27eSBarry Smith 
9691e07b27eSBarry Smith /*@
9701e07b27eSBarry Smith  PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used.
9711e07b27eSBarry Smith 
9721e07b27eSBarry Smith  Not Collective
9731e07b27eSBarry Smith 
9741e07b27eSBarry Smith  Input Parameter:
9751e07b27eSBarry Smith .  pc - the preconditioner context
9761e07b27eSBarry Smith 
9771e07b27eSBarry Smith  Output Parameter:
9781e07b27eSBarry Smith .  v - the flag
9791e07b27eSBarry Smith 
9801e07b27eSBarry Smith  Level: advanced
9811e07b27eSBarry Smith 
9821e07b27eSBarry Smith .keywords: PC, telescoping solve
9831e07b27eSBarry Smith @*/
9841e07b27eSBarry Smith PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v)
9851e07b27eSBarry Smith {
986bd49479cSSatish Balay   PetscErrorCode ierr;
987bd49479cSSatish Balay   PetscFunctionBegin;
988163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr);
989bd49479cSSatish Balay   PetscFunctionReturn(0);
9901e07b27eSBarry Smith }
9911e07b27eSBarry Smith 
9921e07b27eSBarry Smith /*@
9931e07b27eSBarry Smith  PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC.
9941e07b27eSBarry Smith 
9951e07b27eSBarry Smith  Not Collective
9961e07b27eSBarry Smith 
9971e07b27eSBarry Smith  Input Parameter:
9981e07b27eSBarry Smith .  pc - the preconditioner context
9991e07b27eSBarry Smith 
10001e07b27eSBarry Smith  Output Parameter:
10011e07b27eSBarry Smith .  v - Use PETSC_TRUE to ignore any DM
10021e07b27eSBarry Smith 
10031e07b27eSBarry Smith  Level: advanced
10041e07b27eSBarry Smith 
10051e07b27eSBarry Smith .keywords: PC, telescoping solve
10061e07b27eSBarry Smith @*/
1007bfd6bcc6SSatish Balay PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v)
10081e07b27eSBarry Smith {
1009bd49479cSSatish Balay   PetscErrorCode ierr;
1010bd49479cSSatish Balay   PetscFunctionBegin;
1011bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr);
1012bd49479cSSatish Balay   PetscFunctionReturn(0);
10131e07b27eSBarry Smith }
10141e07b27eSBarry Smith 
10151e07b27eSBarry Smith /*@
10168d9f7141SDave May  PCTelescopeGetUseCoarseDM - Get the flag indicating if the coarse DM attached to DM associated with the PC will be used.
10178d9f7141SDave May 
10188d9f7141SDave May  Not Collective
10198d9f7141SDave May 
10208d9f7141SDave May  Input Parameter:
10218d9f7141SDave May .  pc - the preconditioner context
10228d9f7141SDave May 
10238d9f7141SDave May  Output Parameter:
10248d9f7141SDave May .  v - the flag
10258d9f7141SDave May 
10268d9f7141SDave May  Level: advanced
10278d9f7141SDave May 
10288d9f7141SDave May .keywords: PC, telescoping solve
10298d9f7141SDave May @*/
10308d9f7141SDave May PetscErrorCode PCTelescopeGetUseCoarseDM(PC pc,PetscBool *v)
10318d9f7141SDave May {
10328d9f7141SDave May   PetscErrorCode ierr;
10338d9f7141SDave May   PetscFunctionBegin;
10348d9f7141SDave May   ierr = PetscUseMethod(pc,"PCTelescopeGetUseCoarseDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr);
10358d9f7141SDave May   PetscFunctionReturn(0);
10368d9f7141SDave May }
10378d9f7141SDave May 
10388d9f7141SDave May /*@
10398d9f7141SDave May  PCTelescopeSetUseCoarseDM - Set a flag to query the DM attached to the PC if it also has a coarse DM
10408d9f7141SDave May 
10418d9f7141SDave May  Not Collective
10428d9f7141SDave May 
10438d9f7141SDave May  Input Parameter:
10448d9f7141SDave May .  pc - the preconditioner context
10458d9f7141SDave May 
10468d9f7141SDave May  Output Parameter:
10478d9f7141SDave May .  v - Use PETSC_TRUE to ignore any DM
10488d9f7141SDave May 
10498d9f7141SDave May  Notes:
10508d9f7141SDave May  When you have specified to use a coarse DM, the communicator used to create the sub-KSP within PCTelescope
10518d9f7141SDave May  will be that of the coarse DM. Hence the flags -pc_telescope_reduction_factor and
10528d9f7141SDave May  -pc_telescope_subcomm_type will no longer have any meaning.
10538d9f7141SDave May  It is required that the communicator associated with the parent (fine) and the coarse DM are of different sizes.
10548d9f7141SDave May  An error will occur of the size of the communicator associated with the coarse DM
10558d9f7141SDave May  is the same as that of the parent DM.
10568d9f7141SDave May  Furthermore, it is required that the communicator on the coarse DM is a sub-communicator of the parent.
10578d9f7141SDave May  This will be checked at the time the preconditioner is setup and an error will occur if
10588d9f7141SDave May  the coarse DM does not define a sub-communicator of that used by the parent DM.
10598d9f7141SDave May 
10608d9f7141SDave May  The particular Telescope setup invoked when using a coarse DM is agnostic with respect to the type of
10618d9f7141SDave May  the DM used (e.g. it supports DMSHELL, DMPLEX, etc).
10628d9f7141SDave May 
10638d9f7141SDave May  Support is currently only provided for the case when you are using KSPSetComputeOperators()
10648d9f7141SDave May 
10658d9f7141SDave 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.
10668d9f7141SDave May  In the user code, this is achieved via
10678d9f7141SDave May .vb
10688d9f7141SDave May    {
10698d9f7141SDave May      DM dm_fine;
10708d9f7141SDave May      PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeFieldScatter",your_field_scatter_method);
10718d9f7141SDave May    }
10728d9f7141SDave May .ve
10738d9f7141SDave May  The signature of the user provided field scatter method is
10748d9f7141SDave May .vb
10758d9f7141SDave May    PetscErrorCode your_field_scatter_method(DM dm_fine,Vec x_fine,ScatterMode mode,DM dm_coarse,Vec x_coarse);
10768d9f7141SDave May .ve
10778d9f7141SDave May  The user must provide support for both mode = SCATTER_FORWARD and mode = SCATTER_REVERSE.
10788d9f7141SDave May  SCATTER_FORWARD implies the direction of transfer is from the parent (fine) DM to the coarse DM.
10798d9f7141SDave May 
10808d9f7141SDave May  Optionally, the user may also compose a function with the parent DM to facilitate the transfer
10818d9f7141SDave May  of state variables between the fine and coarse DMs.
10828d9f7141SDave May  In the context of a finite element discretization, an example state variable might be
10838d9f7141SDave May  values associated with quadrature points within each element.
10848d9f7141SDave May  A user provided state scatter method is composed via
10858d9f7141SDave May .vb
10868d9f7141SDave May    {
10878d9f7141SDave May      DM dm_fine;
10888d9f7141SDave May      PetscObjectCompose((PetscObject)dm_fine,"PCTelescopeStateScatter",your_state_scatter_method);
10898d9f7141SDave May    }
10908d9f7141SDave May .ve
10918d9f7141SDave May  The signature of the user provided state scatter method is
10928d9f7141SDave May .vb
10938d9f7141SDave May    PetscErrorCode your_state_scatter_method(DM dm_fine,ScatterMode mode,DM dm_coarse);
10948d9f7141SDave May .ve
10958d9f7141SDave May  SCATTER_FORWARD implies the direction of transfer is from the fine DM to the coarse DM.
10968d9f7141SDave May  The user is only required to support mode = SCATTER_FORWARD.
10978d9f7141SDave May  No assumption is made about the data type of the state variables.
10988d9f7141SDave May  These must be managed by the user and must be accessible from the DM.
10998d9f7141SDave May 
11008d9f7141SDave May  Care must be taken in defining the user context passed to KSPSetComputeOperators() which is to be
11018d9f7141SDave May  associated with the sub-KSP residing within PCTelescope.
11028d9f7141SDave May  In general, PCTelescope assumes that the context on the fine and coarse DM used with
11038d9f7141SDave May  KSPSetComputeOperators() should be "similar" in type or origin.
11048d9f7141SDave May  Specifically the following rules are used to infer what context on the sub-KSP should be.
11058d9f7141SDave May 
11068d9f7141SDave May  First the contexts from the KSP and the fine and coarse DMs are retrieved.
11078d9f7141SDave May  Note that the special case of a DMSHELL context is queried.
11088d9f7141SDave May 
11098d9f7141SDave May .vb
11108d9f7141SDave May    DMKSPGetComputeOperators(dm_fine,&dmfine_kspfunc,&dmfine_kspctx);
11118d9f7141SDave May    DMGetApplicationContext(dm_fine,&dmfine_appctx);
11128d9f7141SDave May    DMShellGetContext(dm_fine,&dmfine_shellctx);
11138d9f7141SDave May 
11148d9f7141SDave May    DMGetApplicationContext(dm_coarse,&dmcoarse_appctx);
11158d9f7141SDave May    DMShellGetContext(dm_coarse,&dmcoarse_shellctx);
11168d9f7141SDave May .ve
11178d9f7141SDave May 
11188d9f7141SDave May  The following rules are then enforced:
11198d9f7141SDave May 
11208d9f7141SDave May  1. If dmfine_kspctx = NULL, then we provide a NULL pointer as the context for the sub-KSP:
11218d9f7141SDave May  KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,NULL);
11228d9f7141SDave May 
11238d9f7141SDave May  2. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_appctx,
11248d9f7141SDave May  check that dmcoarse_appctx is also non-NULL. If this is true, then:
11258d9f7141SDave May  KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,dmcoarse_appctx);
11268d9f7141SDave May 
11278d9f7141SDave May  3. If dmfine_kspctx != NULL and dmfine_kspctx == dmfine_shellctx,
11288d9f7141SDave May  check that dmcoarse_shellctx is also non-NULL. If this is true, then:
11298d9f7141SDave May  KSPSetComputeOperators(sub_ksp,dmfine_kspfunc,dmcoarse_shellctx);
11308d9f7141SDave May 
11318d9f7141SDave May  If neither of the above three tests passed, then PCTelescope cannot safely determine what
11328d9f7141SDave May  context should be provided to KSPSetComputeOperators() for use with the sub-KSP.
11338d9f7141SDave May  In this case, an additional mechanism is provided via a composed function which will return
11348d9f7141SDave May  the actual context to be used. To use this feature you must compose the "getter" function
11358d9f7141SDave May  with the coarse DM, e.g.
11368d9f7141SDave May .vb
11378d9f7141SDave May    {
11388d9f7141SDave May      DM dm_coarse;
11398d9f7141SDave May      PetscObjectCompose((PetscObject)dm_coarse,"PCTelescopeGetCoarseDMKSPContext",your_coarse_context_getter);
11408d9f7141SDave May    }
11418d9f7141SDave May .ve
11428d9f7141SDave May  The signature of the user provided method is
11438d9f7141SDave May .vb
11448d9f7141SDave May    PetscErrorCode your_coarse_context_getter(DM dm_coarse,void **your_kspcontext);
11458d9f7141SDave May .ve
11468d9f7141SDave May 
11478d9f7141SDave May  Level: advanced
11488d9f7141SDave May 
11498d9f7141SDave May .keywords: PC, telescoping solve
11508d9f7141SDave May @*/
11518d9f7141SDave May PetscErrorCode PCTelescopeSetUseCoarseDM(PC pc,PetscBool v)
11528d9f7141SDave May {
11538d9f7141SDave May   PetscErrorCode ierr;
11548d9f7141SDave May   PetscFunctionBegin;
11558d9f7141SDave May   ierr = PetscTryMethod(pc,"PCTelescopeSetUseCoarseDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr);
11568d9f7141SDave May   PetscFunctionReturn(0);
11578d9f7141SDave May }
11588d9f7141SDave May 
11598d9f7141SDave May /*@
11600ae7c45bSDave May  PCTelescopeGetIgnoreKSPComputeOperators - Get the flag indicating if KSPComputeOperators will be used.
11610ae7c45bSDave May 
11620ae7c45bSDave May  Not Collective
11630ae7c45bSDave May 
11640ae7c45bSDave May  Input Parameter:
11650ae7c45bSDave May .  pc - the preconditioner context
11660ae7c45bSDave May 
11670ae7c45bSDave May  Output Parameter:
11680ae7c45bSDave May .  v - the flag
11690ae7c45bSDave May 
11700ae7c45bSDave May  Level: advanced
11710ae7c45bSDave May 
11720ae7c45bSDave May .keywords: PC, telescoping solve
11730ae7c45bSDave May @*/
11740ae7c45bSDave May PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC pc,PetscBool *v)
11750ae7c45bSDave May {
11760ae7c45bSDave May   PetscErrorCode ierr;
11770ae7c45bSDave May   PetscFunctionBegin;
1178163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr);
11790ae7c45bSDave May   PetscFunctionReturn(0);
11800ae7c45bSDave May }
11810ae7c45bSDave May 
11820ae7c45bSDave May /*@
11830ae7c45bSDave May  PCTelescopeSetIgnoreKSPComputeOperators - Set a flag to ignore KSPComputeOperators.
11840ae7c45bSDave May 
11850ae7c45bSDave May  Not Collective
11860ae7c45bSDave May 
11870ae7c45bSDave May  Input Parameter:
11880ae7c45bSDave May .  pc - the preconditioner context
11890ae7c45bSDave May 
11900ae7c45bSDave May  Output Parameter:
1191a954d8f4SDave May .  v - Use PETSC_TRUE to ignore the method (if defined) set via KSPSetComputeOperators on pc
11920ae7c45bSDave May 
11930ae7c45bSDave May  Level: advanced
11940ae7c45bSDave May 
11950ae7c45bSDave May .keywords: PC, telescoping solve
11960ae7c45bSDave May @*/
11970ae7c45bSDave May PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC pc,PetscBool v)
11980ae7c45bSDave May {
11990ae7c45bSDave May   PetscErrorCode ierr;
12000ae7c45bSDave May   PetscFunctionBegin;
12010ae7c45bSDave May   ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr);
12020ae7c45bSDave May   PetscFunctionReturn(0);
12030ae7c45bSDave May }
12040ae7c45bSDave May 
12050ae7c45bSDave May /*@
12061e07b27eSBarry Smith  PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP.
12071e07b27eSBarry Smith 
12081e07b27eSBarry Smith  Not Collective
12091e07b27eSBarry Smith 
12101e07b27eSBarry Smith  Input Parameter:
12111e07b27eSBarry Smith .  pc - the preconditioner context
12121e07b27eSBarry Smith 
12131e07b27eSBarry Smith  Output Parameter:
12141e07b27eSBarry Smith .  subdm - The re-partitioned DM
12151e07b27eSBarry Smith 
12161e07b27eSBarry Smith  Level: advanced
12171e07b27eSBarry Smith 
12181e07b27eSBarry Smith .keywords: PC, telescoping solve
12191e07b27eSBarry Smith @*/
12201e07b27eSBarry Smith PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm)
12211e07b27eSBarry Smith {
1222bd49479cSSatish Balay   PetscErrorCode ierr;
1223bd49479cSSatish Balay   PetscFunctionBegin;
1224163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm));CHKERRQ(ierr);
1225bd49479cSSatish Balay   PetscFunctionReturn(0);
12261e07b27eSBarry Smith }
12271e07b27eSBarry Smith 
122848a10b22SPatrick Sanan /*@
122948a10b22SPatrick Sanan  PCTelescopeSetSubcommType - set subcommunicator type (interlaced or contiguous)
123048a10b22SPatrick Sanan 
123148a10b22SPatrick Sanan  Logically Collective
123248a10b22SPatrick Sanan 
123348a10b22SPatrick Sanan  Input Parameter:
12341dae98e4SBarry Smith +  pc - the preconditioner context
12351dae98e4SBarry Smith -  subcommtype - the subcommunicator type (see PetscSubcommType)
123648a10b22SPatrick Sanan 
123748a10b22SPatrick Sanan  Level: advanced
123848a10b22SPatrick Sanan 
123948a10b22SPatrick Sanan .keywords: PC, telescoping solve
124048a10b22SPatrick Sanan 
124148a10b22SPatrick Sanan .seealso: PetscSubcommType, PetscSubcomm, PCTELESCOPE
124248a10b22SPatrick Sanan @*/
124348a10b22SPatrick Sanan PetscErrorCode PCTelescopeSetSubcommType(PC pc, PetscSubcommType subcommtype)
124448a10b22SPatrick Sanan {
124548a10b22SPatrick Sanan   PetscErrorCode ierr;
124648a10b22SPatrick Sanan   PetscFunctionBegin;
124748a10b22SPatrick Sanan   ierr = PetscTryMethod(pc,"PCTelescopeSetSubcommType_C",(PC,PetscSubcommType),(pc,subcommtype));CHKERRQ(ierr);
124848a10b22SPatrick Sanan   PetscFunctionReturn(0);
124948a10b22SPatrick Sanan }
125048a10b22SPatrick Sanan 
125148a10b22SPatrick Sanan /*@
125248a10b22SPatrick Sanan  PCTelescopeGetSubcommType - Get the subcommunicator type (interlaced or contiguous)
125348a10b22SPatrick Sanan 
125448a10b22SPatrick Sanan  Not Collective
125548a10b22SPatrick Sanan 
125648a10b22SPatrick Sanan  Input Parameter:
125748a10b22SPatrick Sanan .  pc - the preconditioner context
125848a10b22SPatrick Sanan 
125948a10b22SPatrick Sanan  Output Parameter:
126048a10b22SPatrick Sanan .  subcommtype - the subcommunicator type (see PetscSubcommType)
126148a10b22SPatrick Sanan 
126248a10b22SPatrick Sanan  Level: advanced
126348a10b22SPatrick Sanan 
126448a10b22SPatrick Sanan .keywords: PC, telescoping solve
126548a10b22SPatrick Sanan 
12661dae98e4SBarry Smith .seealso: PetscSubcomm, PetscSubcommType, PCTELESCOPE
126748a10b22SPatrick Sanan @*/
12681dae98e4SBarry Smith PetscErrorCode PCTelescopeGetSubcommType(PC pc, PetscSubcommType *subcommtype)
126948a10b22SPatrick Sanan {
127048a10b22SPatrick Sanan   PetscErrorCode ierr;
127148a10b22SPatrick Sanan   PetscFunctionBegin;
127248a10b22SPatrick Sanan   ierr = PetscUseMethod(pc,"PCTelescopeGetSubcommType_C",(PC,PetscSubcommType*),(pc,subcommtype));CHKERRQ(ierr);
127348a10b22SPatrick Sanan   PetscFunctionReturn(0);
127448a10b22SPatrick Sanan }
127548a10b22SPatrick Sanan 
12761e07b27eSBarry Smith /* -------------------------------------------------------------------------------------*/
12771e07b27eSBarry Smith /*MC
127800fea0ebSDave May    PCTELESCOPE - Runs a KSP solver on a sub-communicator. MPI ranks not in the sub-communicator are idle during the solve.
12791e07b27eSBarry Smith 
12801e07b27eSBarry Smith    Options Database:
128100fea0ebSDave 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.
128200fea0ebSDave May .  -pc_telescope_ignore_dm  - flag to indicate whether an attached DM should be ignored.
128300fea0ebSDave May .  -pc_telescope_subcomm_type <interlaced,contiguous> - defines the selection of MPI ranks on the sub-communicator. see PetscSubcomm for more information.
128400fea0ebSDave May .  -pc_telescope_ignore_kspcomputeoperators - flag to indicate whether KSPSetComputeOperators should be used on the sub-KSP.
128500fea0ebSDave May -  -pc_telescope_use_coarse_dm - flag to indicate whether the coarse DM should be used to define the sub-communicator.
12861e07b27eSBarry Smith 
12871e07b27eSBarry Smith    Level: advanced
12881e07b27eSBarry Smith 
12891e07b27eSBarry Smith    Notes:
129000fea0ebSDave May    Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation
129100fea0ebSDave May    creates a child sub-communicator (c') containing fewer MPI ranks than the original parent preconditioner (PC).
12926fc41876SBarry Smith    The preconditioner is deemed telescopic as it only calls KSPSolve() on a single
12938439623fSPatrick Sanan    sub-communicator, in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators.
129400fea0ebSDave May    This means there will be MPI ranks which will be idle during the application of this preconditioner.
129500fea0ebSDave May    Additionally, in comparison with PCREDUNDANT, PCTELESCOPE can utilize an attached DM.
12966fc41876SBarry Smith 
129700fea0ebSDave May    The default type of the sub KSP (the KSP defined on c') is PREONLY.
129800fea0ebSDave May 
129900fea0ebSDave May    There are three setup mechanisms for PCTelescope. Features support by each type are described below.
130000fea0ebSDave May    In the following, we will refer to the operators B and B', these are the Bmat provided to the KSP on the
130100fea0ebSDave May    communicators c and c' respectively.
130200fea0ebSDave May 
130300fea0ebSDave May    [1] Default setup
130400fea0ebSDave May    The sub-communicator c' is created via PetscSubcommCreate().
130500fea0ebSDave May    Explicitly defined nullspace and near nullspace vectors will be propogated from B to B'.
130600fea0ebSDave May    Currently there is no support define nullspaces via a user supplied method (e.g. as passed to MatNullSpaceSetFunction()).
130700fea0ebSDave May    No support is provided for KSPSetComputeOperators().
130800fea0ebSDave May    Currently there is no support for the flag -pc_use_amat.
130900fea0ebSDave May 
131000fea0ebSDave May    [2] DM aware setup
131100fea0ebSDave May    If a DM is attached to the PC, it is re-partitioned on the sub-communicator c'.
131200fea0ebSDave May    c' is created via PetscSubcommCreate().
13131e07b27eSBarry Smith    Both the Bmat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned DM.
13141e07b27eSBarry Smith    Currently only support for re-partitioning a DMDA is provided.
131500fea0ebSDave 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').
131600fea0ebSDave May    Currently there is no support define nullspaces via a user supplied method (e.g. as passed to MatNullSpaceSetFunction()).
131700fea0ebSDave May    Support is provided for KSPSetComputeOperators(). The user provided function and context is propagated to the sub KSP.
131800fea0ebSDave May    This is fragile since the user must ensure that their user context is valid for use on c'.
131900fea0ebSDave May    Currently there is no support for the flag -pc_use_amat.
13201e07b27eSBarry Smith 
132100fea0ebSDave May    [3] Coarse DM setup
132200fea0ebSDave May    If a DM (dmfine) is attached to the PC, dmfine is queried for a "coarse" DM (call this dmcoarse) via DMGetCoarseDM().
132300fea0ebSDave May    PCTELESCOPE will interpret the coarse DM as being defined on a sub-communicator of c.
132400fea0ebSDave May    The communicator associated with dmcoarse will define the c' to be used within PCTELESCOPE.
132500fea0ebSDave May    PCTELESCOPE will check that c' is in fact a sub-communicator of c. If it is not, an error will be reported.
132600fea0ebSDave May    The intention of this setup type is that PCTELESCOPE will use an existing (e.g. user defined) communicator hierarchy, say as would be
132700fea0ebSDave May    available with using multi-grid on unstructured meshes.
132800fea0ebSDave May    This setup will not use the command line options -pc_telescope_reduction_factor or -pc_telescope_subcomm_type.
132900fea0ebSDave 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'.
133000fea0ebSDave May    Currently there is no support define nullspaces via a user supplied method (e.g. as passed to MatNullSpaceSetFunction()).
133100fea0ebSDave May    There is no general method to permute field orderings, hence only KSPSetComputeOperators() is supported.
133200fea0ebSDave May    The user must use PetscObjectComposeFunction() with dmfine to define the method to scatter fields from dmfine to dmcoarse.
133300fea0ebSDave 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().
133400fea0ebSDave May    Currently there is no support for the flag -pc_use_amat.
133500fea0ebSDave May    This setup can be invoked by the option -pc_telescope_use_coarse_dm or by calling PCTelescopeSetUseCoarseDM(pc,PETSC_TRUE);
133600fea0ebSDave May    Further information about the user-provided methods required by this setup type are described here PCTelescopeSetUseCoarseDM().
13376fc41876SBarry Smith 
13386fc41876SBarry Smith    Developer Notes:
13396fc41876SBarry Smith    During PCSetup, the B operator is scattered onto c'.
13406fc41876SBarry Smith    Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c').
13418439623fSPatrick Sanan    Then, KSPSolve() is executed on the c' communicator.
13426fc41876SBarry Smith 
13436fc41876SBarry Smith    The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED
1344a04a6428SPatrick 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.
13456fc41876SBarry Smith 
1346005d9f20SPatrick Sanan    The telescoping preconditioner is aware of nullspaces and near nullspaces which are attached to the B operator.
13478439623fSPatrick Sanan    In the case where B has a (near) nullspace attached, the (near) nullspace vectors are extracted from B and mapped into
1348005d9f20SPatrick Sanan    a new (near) nullspace, defined on the sub-communicator, which is attached to B' (the B operator which was scattered to c')
13496fc41876SBarry Smith 
135000fea0ebSDave May    The telescoping preconditioner can re-partition an attached DM if it is a DMDA (2D or 3D -
135100fea0ebSDave May    support for 1D DMDAs is not provided). If a DMDA is found, a topolgically equivalent DMDA is created on c'
135200fea0ebSDave May    and this new DM is attached the sub KSP. The design of telescope is such that it should be possible to extend support
13538439623fSPatrick Sanan    for re-partitioning other to DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs.
135400fea0ebSDave May    Alternatively, user-provided re-partitioned DMs can be used via -pc_telescope_use_coarse_dm.
13556fc41876SBarry Smith 
135600fea0ebSDave May    With the default setup mode, B' is defined by fusing rows (in order) associated with MPI ranks common to c and c'.
13576fc41876SBarry Smith 
13588439623fSPatrick Sanan    When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permutation of B
13597dae84e0SHong Zhang    into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatCreateSubMatrices(), (iii) fusing the
13606fc41876SBarry Smith    locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat()
13616fc41876SBarry Smith 
13628439623fSPatrick Sanan    Limitations/improvements include the following.
13638439623fSPatrick Sanan    VecPlaceArray() could be used within PCApply() to improve efficiency and reduce memory usage.
136400fea0ebSDave May    A unified mechanism to query for user contexts as required by KSPSetComputeOperators() and MatNullSpaceSetFunction().
13656fc41876SBarry Smith 
13666fc41876SBarry Smith    The symmetric permutation used when a DMDA is encountered is performed via explicitly assmbleming a permutation matrix P,
13678439623fSPatrick 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
13688439623fSPatrick Sanan    VecPermute() does not supported for the use case required here. By computing P, one can permute both the operator and RHS in a
13696fc41876SBarry Smith    consistent manner.
13706fc41876SBarry Smith 
137100fea0ebSDave May    Mapping of vectors (default setup mode) is performed in the following way.
137200fea0ebSDave 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.
13738439623fSPatrick Sanan    Using the interlaced creation routine, the ranks in c with color = 0 will be rank 0 and 2.
137400fea0ebSDave May    We perform the scatter to the sub-communicator in the following way.
1375514db83aSDave May    [1] Given a vector x defined on communicator c
13766fc41876SBarry Smith 
1377514db83aSDave May .vb
1378514db83aSDave May    rank(c)  local values of x
1379514db83aSDave May    ------- ----------------------------------------
1380514db83aSDave May         0   [  0.0,  1.0,  2.0,  3.0,  4.0,  5.0 ]
1381514db83aSDave May         1   [  6.0,  7.0,  8.0,  9.0, 10.0, 11.0 ]
1382514db83aSDave May         2   [ 12.0, 13.0, 14.0, 15.0, 16.0, 17.0 ]
1383514db83aSDave May         3   [ 18.0, 19.0, 20.0, 21.0, 22.0, 23.0 ]
1384514db83aSDave May .ve
13856fc41876SBarry Smith 
1386514db83aSDave May    scatter into xtmp defined also on comm c, so that we have the following values
13876fc41876SBarry Smith 
1388514db83aSDave May .vb
1389514db83aSDave May    rank(c)  local values of xtmp
1390514db83aSDave May    ------- ----------------------------------------------------------------------------
1391514db83aSDave 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 ]
1392514db83aSDave May         1   [ ]
1393514db83aSDave 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 ]
1394514db83aSDave May         3   [ ]
1395514db83aSDave May .ve
13966fc41876SBarry Smith 
13976fc41876SBarry Smith    The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') have no values
13986fc41876SBarry Smith 
13996fc41876SBarry Smith 
1400514db83aSDave 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'.
14016fc41876SBarry Smith    Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0.
14026fc41876SBarry Smith 
1403514db83aSDave May .vb
1404514db83aSDave May    rank(c')  local values of xred
1405514db83aSDave May    -------- ----------------------------------------------------------------------------
1406514db83aSDave 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 ]
1407514db83aSDave 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 ]
1408514db83aSDave May .ve
14091e07b27eSBarry Smith 
14101e07b27eSBarry Smith   Contributed by Dave May
14111e07b27eSBarry Smith 
1412bf00f589SPatrick Sanan   Reference:
1413bf00f589SPatrick 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
1414bf00f589SPatrick Sanan 
14156fc41876SBarry Smith .seealso:  PCTelescopeGetKSP(), PCTelescopeGetDM(), PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(), PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT
14161e07b27eSBarry Smith M*/
14171e07b27eSBarry Smith PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc)
14181e07b27eSBarry Smith {
14191e07b27eSBarry Smith   PetscErrorCode       ierr;
14201e07b27eSBarry Smith   struct _PC_Telescope *sred;
14211e07b27eSBarry Smith 
14221e07b27eSBarry Smith   PetscFunctionBegin;
14231e07b27eSBarry Smith   ierr = PetscNewLog(pc,&sred);CHKERRQ(ierr);
14242a22aa42SDave May   sred->psubcomm       = NULL;
142548a10b22SPatrick Sanan   sred->subcommtype    = PETSC_SUBCOMM_INTERLACED;
14262a22aa42SDave May   sred->subcomm        = MPI_COMM_NULL;
14271e07b27eSBarry Smith   sred->redfactor      = 1;
14281e07b27eSBarry Smith   sred->ignore_dm      = PETSC_FALSE;
14297c5279cbSDave May   sred->ignore_kspcomputeoperators = PETSC_FALSE;
14308d9f7141SDave May   sred->use_coarse_dm  = PETSC_FALSE;
14311e07b27eSBarry Smith   pc->data             = (void*)sred;
14321e07b27eSBarry Smith 
14331e07b27eSBarry Smith   pc->ops->apply           = PCApply_Telescope;
14341e07b27eSBarry Smith   pc->ops->applytranspose  = NULL;
1435f650675bSDave May   pc->ops->applyrichardson = PCApplyRichardson_Telescope;
14361e07b27eSBarry Smith   pc->ops->setup           = PCSetUp_Telescope;
14371e07b27eSBarry Smith   pc->ops->destroy         = PCDestroy_Telescope;
14381e07b27eSBarry Smith   pc->ops->reset           = PCReset_Telescope;
14391e07b27eSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Telescope;
14401e07b27eSBarry Smith   pc->ops->view            = PCView_Telescope;
14411e07b27eSBarry Smith 
14421e07b27eSBarry Smith   sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
14431e07b27eSBarry Smith   sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
14441e07b27eSBarry Smith   sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
14451e07b27eSBarry Smith   sred->pctelescope_reset_type              = NULL;
14461e07b27eSBarry Smith 
14471e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetKSP_C",PCTelescopeGetKSP_Telescope);CHKERRQ(ierr);
144848a10b22SPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetSubcommType_C",PCTelescopeGetSubcommType_Telescope);CHKERRQ(ierr);
144948a10b22SPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetSubcommType_C",PCTelescopeSetSubcommType_Telescope);CHKERRQ(ierr);
14501e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetReductionFactor_C",PCTelescopeGetReductionFactor_Telescope);CHKERRQ(ierr);
14511e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetReductionFactor_C",PCTelescopeSetReductionFactor_Telescope);CHKERRQ(ierr);
14521e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreDM_C",PCTelescopeGetIgnoreDM_Telescope);CHKERRQ(ierr);
14531e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreDM_C",PCTelescopeSetIgnoreDM_Telescope);CHKERRQ(ierr);
14540ae7c45bSDave May   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",PCTelescopeGetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr);
14550ae7c45bSDave May   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",PCTelescopeSetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr);
14561e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetDM_C",PCTelescopeGetDM_Telescope);CHKERRQ(ierr);
14578d9f7141SDave May   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetUseCoarseDM_C",PCTelescopeGetUseCoarseDM_Telescope);CHKERRQ(ierr);
14588d9f7141SDave May   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetUseCoarseDM_C",PCTelescopeSetUseCoarseDM_Telescope);CHKERRQ(ierr);
14591e07b27eSBarry Smith   PetscFunctionReturn(0);
14601e07b27eSBarry Smith }
1461