xref: /petsc/src/ksp/pc/impls/telescope/telescope.c (revision 2a22aa427fba2a1731631377d9475e4b679dc0e4)
11e07b27eSBarry Smith 
2120bdd93SDave May #include <petsc/private/matimpl.h>
36fc41876SBarry Smith #include <petsc/private/pcimpl.h>
41e07b27eSBarry Smith #include <petscksp.h> /*I "petscksp.h" I*/
51e07b27eSBarry Smith #include <petscdm.h> /*I "petscdm.h" I*/
6575a0592SBarry Smith #include "../src/ksp/pc/impls/telescope/telescope.h"
71e07b27eSBarry Smith 
8bf00f589SPatrick Sanan static PetscBool  cited = PETSC_FALSE;
9bf00f589SPatrick Sanan static const char citation[] =
10bf00f589SPatrick Sanan "@inproceedings{MaySananRuppKnepleySmith2016,\n"
11bf00f589SPatrick Sanan "  title     = {Extreme-Scale Multigrid Components within PETSc},\n"
12bf00f589SPatrick Sanan "  author    = {Dave A. May and Patrick Sanan and Karl Rupp and Matthew G. Knepley and Barry F. Smith},\n"
13bf00f589SPatrick Sanan "  booktitle = {Proceedings of the Platform for Advanced Scientific Computing Conference},\n"
14bf00f589SPatrick Sanan "  series    = {PASC '16},\n"
15bf00f589SPatrick Sanan "  isbn      = {978-1-4503-4126-4},\n"
16bf00f589SPatrick Sanan "  location  = {Lausanne, Switzerland},\n"
17bf00f589SPatrick Sanan "  pages     = {5:1--5:12},\n"
18bf00f589SPatrick Sanan "  articleno = {5},\n"
19bf00f589SPatrick Sanan "  numpages  = {12},\n"
20bf00f589SPatrick Sanan "  url       = {http://doi.acm.org/10.1145/2929908.2929913},\n"
21bf00f589SPatrick Sanan "  doi       = {10.1145/2929908.2929913},\n"
22bf00f589SPatrick Sanan "  acmid     = {2929913},\n"
23bf00f589SPatrick Sanan "  publisher = {ACM},\n"
24bf00f589SPatrick Sanan "  address   = {New York, NY, USA},\n"
25bf00f589SPatrick Sanan "  keywords  = {GPU, HPC, agglomeration, coarse-level solver, multigrid, parallel computing, preconditioning},\n"
26bf00f589SPatrick Sanan "  year      = {2016}\n"
27bf00f589SPatrick Sanan "}\n";
28bf00f589SPatrick Sanan 
291e07b27eSBarry Smith /*
301e07b27eSBarry Smith  PCTelescopeSetUp_default()
311e07b27eSBarry Smith  PCTelescopeMatCreate_default()
321e07b27eSBarry Smith 
331e07b27eSBarry Smith  default
341e07b27eSBarry Smith 
351e07b27eSBarry Smith  // scatter in
361e07b27eSBarry Smith  x(comm) -> xtmp(comm)
371e07b27eSBarry Smith 
381e07b27eSBarry Smith  xred(subcomm) <- xtmp
391e07b27eSBarry Smith  yred(subcomm)
401e07b27eSBarry Smith 
411e07b27eSBarry Smith  yred(subcomm) --> xtmp
421e07b27eSBarry Smith 
431e07b27eSBarry Smith  // scatter out
441e07b27eSBarry Smith  xtmp(comm) -> y(comm)
451e07b27eSBarry Smith */
461e07b27eSBarry Smith 
47*2a22aa42SDave May PetscBool PetscSubComm_isActiveRank(PetscSubcomm scomm)
481e07b27eSBarry Smith {
491e07b27eSBarry Smith   if (scomm->color == 0) { return PETSC_TRUE; }
501e07b27eSBarry Smith   else { return PETSC_FALSE; }
511e07b27eSBarry Smith }
521e07b27eSBarry Smith 
53*2a22aa42SDave May PetscBool isActiveRank(PC_Telescope sred)
54*2a22aa42SDave May {
55*2a22aa42SDave May   if (sred->psubcomm) {
56*2a22aa42SDave May     return(PetscSubComm_isActiveRank(sred->psubcomm));
57*2a22aa42SDave May   } else {
58*2a22aa42SDave May     if (sred->subcomm != MPI_COMM_NULL) { return PETSC_TRUE; }
59*2a22aa42SDave May     else { return PETSC_FALSE; }
60*2a22aa42SDave May   }
61*2a22aa42SDave May }
62*2a22aa42SDave May 
631e07b27eSBarry Smith DM private_PCTelescopeGetSubDM(PC_Telescope sred)
641e07b27eSBarry Smith {
65c6a0d831SBarry Smith   DM subdm = NULL;
661e07b27eSBarry Smith 
67*2a22aa42SDave May   if (!isActiveRank(sred)) { subdm = NULL; }
681e07b27eSBarry Smith   else {
691e07b27eSBarry Smith     switch (sred->sr_type) {
701e07b27eSBarry Smith     case TELESCOPE_DEFAULT: subdm = NULL;
711e07b27eSBarry Smith       break;
721e07b27eSBarry Smith     case TELESCOPE_DMDA:    subdm = ((PC_Telescope_DMDACtx*)sred->dm_ctx)->dmrepart;
731e07b27eSBarry Smith       break;
741e07b27eSBarry Smith     case TELESCOPE_DMPLEX:  subdm = NULL;
751e07b27eSBarry Smith       break;
761e07b27eSBarry Smith     }
771e07b27eSBarry Smith   }
781e07b27eSBarry Smith   return(subdm);
791e07b27eSBarry Smith }
801e07b27eSBarry Smith 
811e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_default(PC pc,PC_Telescope sred)
821e07b27eSBarry Smith {
831e07b27eSBarry Smith   PetscErrorCode ierr;
841e07b27eSBarry Smith   PetscInt       m,M,bs,st,ed;
851e07b27eSBarry Smith   Vec            x,xred,yred,xtmp;
861e07b27eSBarry Smith   Mat            B;
871e07b27eSBarry Smith   MPI_Comm       comm,subcomm;
881e07b27eSBarry Smith   VecScatter     scatter;
891e07b27eSBarry Smith   IS             isin;
901e07b27eSBarry Smith 
911e07b27eSBarry Smith   PetscFunctionBegin;
921e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: setup (default)\n");CHKERRQ(ierr);
931e07b27eSBarry Smith   comm = PetscSubcommParent(sred->psubcomm);
941e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
951e07b27eSBarry Smith 
961e07b27eSBarry Smith   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
971e07b27eSBarry Smith   ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr);
981e07b27eSBarry Smith   ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr);
991e07b27eSBarry Smith   ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr);
1001e07b27eSBarry Smith 
1011e07b27eSBarry Smith   xred = NULL;
1023ac26c5eSBarry Smith   m    = 0;
103*2a22aa42SDave May   if (isActiveRank(sred)) {
1041e07b27eSBarry Smith     ierr = VecCreate(subcomm,&xred);CHKERRQ(ierr);
1051e07b27eSBarry Smith     ierr = VecSetSizes(xred,PETSC_DECIDE,M);CHKERRQ(ierr);
1061e07b27eSBarry Smith     ierr = VecSetBlockSize(xred,bs);CHKERRQ(ierr);
1071e07b27eSBarry Smith     ierr = VecSetFromOptions(xred);CHKERRQ(ierr);
108ca43db0aSBarry Smith     ierr = VecGetLocalSize(xred,&m);CHKERRQ(ierr);
1091e07b27eSBarry Smith   }
1101e07b27eSBarry Smith 
1111e07b27eSBarry Smith   yred = NULL;
112*2a22aa42SDave May   if (isActiveRank(sred)) {
1131e07b27eSBarry Smith     ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr);
1141e07b27eSBarry Smith   }
1151e07b27eSBarry Smith 
1161e07b27eSBarry Smith   ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr);
1171e07b27eSBarry Smith   ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr);
1181e07b27eSBarry Smith   ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr);
1191e07b27eSBarry Smith   ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr);
1201e07b27eSBarry Smith 
121*2a22aa42SDave May   if (isActiveRank(sred)) {
1221e07b27eSBarry Smith     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
1231e07b27eSBarry Smith     ierr = ISCreateStride(comm,(ed-st),st,1,&isin);CHKERRQ(ierr);
1241e07b27eSBarry Smith   } else {
1251e07b27eSBarry Smith     ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr);
1263ac26c5eSBarry Smith     ierr = ISCreateStride(comm,0,st,1,&isin);CHKERRQ(ierr);
1271e07b27eSBarry Smith   }
1281e07b27eSBarry Smith   ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr);
1291e07b27eSBarry Smith 
13035928de7SBarry Smith   ierr = VecScatterCreateWithData(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr);
1311e07b27eSBarry Smith 
1321e07b27eSBarry Smith   sred->isin    = isin;
1331e07b27eSBarry Smith   sred->scatter = scatter;
1341e07b27eSBarry Smith   sred->xred    = xred;
1351e07b27eSBarry Smith   sred->yred    = yred;
1361e07b27eSBarry Smith   sred->xtmp    = xtmp;
1371e07b27eSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
1381e07b27eSBarry Smith   PetscFunctionReturn(0);
1391e07b27eSBarry Smith }
1401e07b27eSBarry Smith 
1411e07b27eSBarry Smith PetscErrorCode PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A)
1421e07b27eSBarry Smith {
1431e07b27eSBarry Smith   PetscErrorCode ierr;
1441e07b27eSBarry Smith   MPI_Comm       comm,subcomm;
1451e07b27eSBarry Smith   Mat            Bred,B;
1461e07b27eSBarry Smith   PetscInt       nr,nc;
1471e07b27eSBarry Smith   IS             isrow,iscol;
1481e07b27eSBarry Smith   Mat            Blocal,*_Blocal;
1491e07b27eSBarry Smith 
1501e07b27eSBarry Smith   PetscFunctionBegin;
1511e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (default)\n");CHKERRQ(ierr);
1521e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1531e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
1541e07b27eSBarry Smith   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
1551e07b27eSBarry Smith   ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr);
1561e07b27eSBarry Smith   isrow = sred->isin;
1571e07b27eSBarry Smith   ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr);
1587dae84e0SHong Zhang   ierr = MatCreateSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr);
1591e07b27eSBarry Smith   Blocal = *_Blocal;
1601e07b27eSBarry Smith   ierr = PetscFree(_Blocal);CHKERRQ(ierr);
1611e07b27eSBarry Smith   Bred = NULL;
162*2a22aa42SDave May   if (isActiveRank(sred)) {
1631e07b27eSBarry Smith     PetscInt mm;
1641e07b27eSBarry Smith 
1651e07b27eSBarry Smith     if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; }
1661e07b27eSBarry Smith 
1671e07b27eSBarry Smith     ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr);
1681e07b27eSBarry Smith     ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr);
1691e07b27eSBarry Smith   }
1701e07b27eSBarry Smith   *A = Bred;
1711e07b27eSBarry Smith   ierr = ISDestroy(&iscol);CHKERRQ(ierr);
1721e07b27eSBarry Smith   ierr = MatDestroy(&Blocal);CHKERRQ(ierr);
1731e07b27eSBarry Smith   PetscFunctionReturn(0);
1741e07b27eSBarry Smith }
1751e07b27eSBarry Smith 
176392968a1SPatrick Sanan static PetscErrorCode PCTelescopeSubNullSpaceCreate_Telescope(PC pc,PC_Telescope sred,MatNullSpace nullspace,MatNullSpace *sub_nullspace)
1771e07b27eSBarry Smith {
1781e07b27eSBarry Smith   PetscErrorCode   ierr;
1791e07b27eSBarry Smith   PetscBool        has_const;
1801e07b27eSBarry Smith   const Vec        *vecs;
181c41e779fSDave May   Vec              *sub_vecs = NULL;
182392968a1SPatrick Sanan   PetscInt         i,k,n = 0;
1831e07b27eSBarry Smith   MPI_Comm         subcomm;
1841e07b27eSBarry Smith 
1851e07b27eSBarry Smith   PetscFunctionBegin;
1861e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
1871e07b27eSBarry Smith   ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr);
1881e07b27eSBarry Smith 
189*2a22aa42SDave May   if (isActiveRank(sred)) {
190e3acf2f7SBarry Smith     if (n) {
191e3acf2f7SBarry Smith       ierr = VecDuplicateVecs(sred->xred,n,&sub_vecs);CHKERRQ(ierr);
1921e07b27eSBarry Smith     }
1931e07b27eSBarry Smith   }
1941e07b27eSBarry Smith 
1951e07b27eSBarry Smith   /* copy entries */
1961e07b27eSBarry Smith   for (k=0; k<n; k++) {
1971e07b27eSBarry Smith     const PetscScalar *x_array;
1981e07b27eSBarry Smith     PetscScalar       *LA_sub_vec;
19913c30530SDave May     PetscInt          st,ed;
2001e07b27eSBarry Smith 
2011e07b27eSBarry Smith     /* pull in vector x->xtmp */
2021e07b27eSBarry Smith     ierr = VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2031e07b27eSBarry Smith     ierr = VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
20447856c66SBarry Smith     if (sub_vecs) {
205a04a6428SPatrick Sanan       /* copy vector entries into xred */
2061e07b27eSBarry Smith       ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
207ea2b237eSDave May       if (sub_vecs[k]) {
2081e07b27eSBarry Smith         ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr);
2091e07b27eSBarry Smith         ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
2101e07b27eSBarry Smith         for (i=0; i<ed-st; i++) {
2111e07b27eSBarry Smith           LA_sub_vec[i] = x_array[i];
2121e07b27eSBarry Smith         }
2131e07b27eSBarry Smith         ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
2141e07b27eSBarry Smith       }
2151e07b27eSBarry Smith       ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
2161e07b27eSBarry Smith     }
21747856c66SBarry Smith   }
2181e07b27eSBarry Smith 
219*2a22aa42SDave May   if (isActiveRank(sred)) {
220d8b9d5b7SPatrick Sanan     /* create new (near) nullspace for redundant object */
221392968a1SPatrick Sanan     ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,sub_nullspace);CHKERRQ(ierr);
222392968a1SPatrick Sanan     ierr = VecDestroyVecs(n,&sub_vecs);CHKERRQ(ierr);
223d8b9d5b7SPatrick Sanan     if (nullspace->remove) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope");
224d8b9d5b7SPatrick 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");
225d8b9d5b7SPatrick Sanan   }
226392968a1SPatrick Sanan   PetscFunctionReturn(0);
227392968a1SPatrick Sanan }
228392968a1SPatrick Sanan 
229392968a1SPatrick Sanan static PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat)
230392968a1SPatrick Sanan {
231392968a1SPatrick Sanan   PetscErrorCode   ierr;
232392968a1SPatrick Sanan   Mat              B;
233392968a1SPatrick Sanan 
234392968a1SPatrick Sanan   PetscFunctionBegin;
235392968a1SPatrick Sanan   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
236392968a1SPatrick Sanan 
237392968a1SPatrick Sanan   /* Propagate the nullspace if it exists */
238392968a1SPatrick Sanan   {
239392968a1SPatrick Sanan     MatNullSpace nullspace,sub_nullspace;
240392968a1SPatrick Sanan     ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr);
241392968a1SPatrick Sanan     if (nullspace) {
242392968a1SPatrick Sanan       ierr = PetscInfo(pc,"PCTelescope: generating nullspace (default)\n");CHKERRQ(ierr);
243392968a1SPatrick Sanan       ierr = PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nullspace,&sub_nullspace);CHKERRQ(ierr);
244*2a22aa42SDave May       if (isActiveRank(sred)) {
245392968a1SPatrick Sanan         ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr);
24641ff1ee9SPatrick Sanan         ierr = MatNullSpaceDestroy(&sub_nullspace);CHKERRQ(ierr);
2471e07b27eSBarry Smith       }
248392968a1SPatrick Sanan     }
249392968a1SPatrick Sanan   }
250392968a1SPatrick Sanan 
251392968a1SPatrick Sanan   /* Propagate the near nullspace if it exists */
252392968a1SPatrick Sanan   {
253392968a1SPatrick Sanan     MatNullSpace nearnullspace,sub_nearnullspace;
254392968a1SPatrick Sanan     ierr = MatGetNearNullSpace(B,&nearnullspace);CHKERRQ(ierr);
255392968a1SPatrick Sanan     if (nearnullspace) {
256392968a1SPatrick Sanan       ierr = PetscInfo(pc,"PCTelescope: generating near nullspace (default)\n");CHKERRQ(ierr);
257392968a1SPatrick Sanan       ierr = PCTelescopeSubNullSpaceCreate_Telescope(pc,sred,nearnullspace,&sub_nearnullspace);CHKERRQ(ierr);
258*2a22aa42SDave May       if (isActiveRank(sred)) {
259392968a1SPatrick Sanan         ierr = MatSetNearNullSpace(sub_mat,sub_nearnullspace);CHKERRQ(ierr);
260392968a1SPatrick Sanan         ierr = MatNullSpaceDestroy(&sub_nearnullspace);CHKERRQ(ierr);
261392968a1SPatrick Sanan       }
262392968a1SPatrick Sanan     }
263392968a1SPatrick Sanan   }
2641e07b27eSBarry Smith   PetscFunctionReturn(0);
2651e07b27eSBarry Smith }
2661e07b27eSBarry Smith 
2671e07b27eSBarry Smith static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer)
2681e07b27eSBarry Smith {
2691e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
2701e07b27eSBarry Smith   PetscErrorCode   ierr;
2711e07b27eSBarry Smith   PetscBool        iascii,isstring;
2721e07b27eSBarry Smith   PetscViewer      subviewer;
2731e07b27eSBarry Smith 
2741e07b27eSBarry Smith   PetscFunctionBegin;
2751e07b27eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
2761e07b27eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
2771e07b27eSBarry Smith   if (iascii) {
2781e07b27eSBarry Smith     if (!sred->psubcomm) {
279efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  preconditioner not yet setup\n");CHKERRQ(ierr);
2801e07b27eSBarry Smith     } else {
2811e07b27eSBarry Smith       MPI_Comm    comm,subcomm;
2821e07b27eSBarry Smith       PetscMPIInt comm_size,subcomm_size;
2831e07b27eSBarry Smith       DM          dm,subdm;
2841e07b27eSBarry Smith 
2851e07b27eSBarry Smith       ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
2861e07b27eSBarry Smith       subdm = private_PCTelescopeGetSubDM(sred);
2871e07b27eSBarry Smith       comm = PetscSubcommParent(sred->psubcomm);
2881e07b27eSBarry Smith       subcomm = PetscSubcommChild(sred->psubcomm);
2891e07b27eSBarry Smith       ierr = MPI_Comm_size(comm,&comm_size);CHKERRQ(ierr);
2901e07b27eSBarry Smith       ierr = MPI_Comm_size(subcomm,&subcomm_size);CHKERRQ(ierr);
2911e07b27eSBarry Smith 
292efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  parent comm size reduction factor = %D\n",sred->redfactor);CHKERRQ(ierr);
293efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  comm_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size);CHKERRQ(ierr);
29448a10b22SPatrick Sanan       switch (sred->subcommtype) {
29548a10b22SPatrick Sanan         case PETSC_SUBCOMM_INTERLACED :
296efd4aadfSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"  subcomm type: interlaced\n",sred->subcommtype);CHKERRQ(ierr);
29748a10b22SPatrick Sanan           break;
29848a10b22SPatrick Sanan         case PETSC_SUBCOMM_CONTIGUOUS :
299efd4aadfSBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,"  subcomm type: contiguous\n",sred->subcommtype);CHKERRQ(ierr);
30048a10b22SPatrick Sanan           break;
30148a10b22SPatrick Sanan         default :
30248a10b22SPatrick Sanan           SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"General subcomm type not supported by PCTelescope");
30348a10b22SPatrick Sanan       }
3041e07b27eSBarry Smith       ierr = PetscViewerGetSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
305*2a22aa42SDave May       if (isActiveRank(sred)) {
3061e07b27eSBarry Smith         ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr);
3071e07b27eSBarry Smith 
3081e07b27eSBarry Smith         if (dm && sred->ignore_dm) {
309efd4aadfSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  ignoring DM\n");CHKERRQ(ierr);
3101e07b27eSBarry Smith         }
3117c5279cbSDave May         if (sred->ignore_kspcomputeoperators) {
312efd4aadfSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  ignoring KSPComputeOperators\n");CHKERRQ(ierr);
3137c5279cbSDave May         }
3141e07b27eSBarry Smith         switch (sred->sr_type) {
3151e07b27eSBarry Smith         case TELESCOPE_DEFAULT:
316efd4aadfSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  using default setup\n");CHKERRQ(ierr);
3171e07b27eSBarry Smith           break;
3181e07b27eSBarry Smith         case TELESCOPE_DMDA:
319efd4aadfSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  DMDA detected\n");CHKERRQ(ierr);
3208ef9ca65SPatrick Sanan           ierr = DMView_DA_Short(subdm,subviewer);CHKERRQ(ierr);
3211e07b27eSBarry Smith           break;
3221e07b27eSBarry Smith         case TELESCOPE_DMPLEX:
323efd4aadfSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  DMPLEX detected\n");CHKERRQ(ierr);
3241e07b27eSBarry Smith           break;
3251e07b27eSBarry Smith         }
3261e07b27eSBarry Smith 
3271e07b27eSBarry Smith         ierr = KSPView(sred->ksp,subviewer);CHKERRQ(ierr);
3281e07b27eSBarry Smith         ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr);
3291e07b27eSBarry Smith       }
3301e07b27eSBarry Smith       ierr = PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
3311e07b27eSBarry Smith     }
3321e07b27eSBarry Smith   }
3331e07b27eSBarry Smith   PetscFunctionReturn(0);
3341e07b27eSBarry Smith }
3351e07b27eSBarry Smith 
3361e07b27eSBarry Smith static PetscErrorCode PCSetUp_Telescope(PC pc)
3371e07b27eSBarry Smith {
3381e07b27eSBarry Smith   PC_Telescope      sred = (PC_Telescope)pc->data;
3391e07b27eSBarry Smith   PetscErrorCode    ierr;
340bd49479cSSatish Balay   MPI_Comm          comm,subcomm=0;
3411e07b27eSBarry Smith   PCTelescopeType   sr_type;
3421e07b27eSBarry Smith 
3431e07b27eSBarry Smith   PetscFunctionBegin;
3441e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
3451e07b27eSBarry Smith 
3461e07b27eSBarry Smith   /* subcomm definition */
3471e07b27eSBarry Smith   if (!pc->setupcalled) {
3481e07b27eSBarry Smith     if (!sred->psubcomm) {
3491e07b27eSBarry Smith       ierr = PetscSubcommCreate(comm,&sred->psubcomm);CHKERRQ(ierr);
3501e07b27eSBarry Smith       ierr = PetscSubcommSetNumber(sred->psubcomm,sred->redfactor);CHKERRQ(ierr);
35148a10b22SPatrick Sanan       ierr = PetscSubcommSetType(sred->psubcomm,sred->subcommtype);CHKERRQ(ierr);
3521e07b27eSBarry Smith       ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr);
353*2a22aa42SDave May       sred->subcomm = PetscSubcommChild(sred->psubcomm);
3541e07b27eSBarry Smith     }
3551e07b27eSBarry Smith   }
356*2a22aa42SDave May   subcomm = sred->subcomm;
3571e07b27eSBarry Smith 
3581e07b27eSBarry Smith   /* internal KSP */
3591e07b27eSBarry Smith   if (!pc->setupcalled) {
3601e07b27eSBarry Smith     const char *prefix;
3611e07b27eSBarry Smith 
362*2a22aa42SDave May     if (isActiveRank(sred)) {
3631e07b27eSBarry Smith       ierr = KSPCreate(subcomm,&sred->ksp);CHKERRQ(ierr);
3641e07b27eSBarry Smith       ierr = KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure);CHKERRQ(ierr);
3651e07b27eSBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
3661e07b27eSBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp);CHKERRQ(ierr);
3671e07b27eSBarry Smith       ierr = KSPSetType(sred->ksp,KSPPREONLY);CHKERRQ(ierr);
3681e07b27eSBarry Smith       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
3691e07b27eSBarry Smith       ierr = KSPSetOptionsPrefix(sred->ksp,prefix);CHKERRQ(ierr);
3701e07b27eSBarry Smith       ierr = KSPAppendOptionsPrefix(sred->ksp,"telescope_");CHKERRQ(ierr);
3711e07b27eSBarry Smith     }
3721e07b27eSBarry Smith   }
3731e07b27eSBarry Smith   /* Determine type of setup/update */
3741e07b27eSBarry Smith   if (!pc->setupcalled) {
3751e07b27eSBarry Smith     PetscBool has_dm,same;
3761e07b27eSBarry Smith     DM        dm;
3771e07b27eSBarry Smith 
3781e07b27eSBarry Smith     sr_type = TELESCOPE_DEFAULT;
3791e07b27eSBarry Smith     has_dm = PETSC_FALSE;
3801e07b27eSBarry Smith     ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
3811e07b27eSBarry Smith     if (dm) { has_dm = PETSC_TRUE; }
3821e07b27eSBarry Smith     if (has_dm) {
3831e07b27eSBarry Smith       /* check for dmda */
3841e07b27eSBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&same);CHKERRQ(ierr);
3851e07b27eSBarry Smith       if (same) {
3861e07b27eSBarry Smith         ierr = PetscInfo(pc,"PCTelescope: found DMDA\n");CHKERRQ(ierr);
3871e07b27eSBarry Smith         sr_type = TELESCOPE_DMDA;
3881e07b27eSBarry Smith       }
3891e07b27eSBarry Smith       /* check for dmplex */
3901e07b27eSBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same);CHKERRQ(ierr);
3911e07b27eSBarry Smith       if (same) {
392994fe344SLisandro Dalcin         ierr = PetscInfo(pc,"PCTelescope: found DMPLEX\n");CHKERRQ(ierr);
3931e07b27eSBarry Smith         sr_type = TELESCOPE_DMPLEX;
3941e07b27eSBarry Smith       }
3951e07b27eSBarry Smith     }
3961e07b27eSBarry Smith 
3971e07b27eSBarry Smith     if (sred->ignore_dm) {
398994fe344SLisandro Dalcin       ierr = PetscInfo(pc,"PCTelescope: ignore DM\n");CHKERRQ(ierr);
3991e07b27eSBarry Smith       sr_type = TELESCOPE_DEFAULT;
4001e07b27eSBarry Smith     }
4011e07b27eSBarry Smith     sred->sr_type = sr_type;
4021e07b27eSBarry Smith   } else {
4031e07b27eSBarry Smith     sr_type = sred->sr_type;
4041e07b27eSBarry Smith   }
4051e07b27eSBarry Smith 
406d8b9d5b7SPatrick Sanan   /* set function pointers for repartition setup, matrix creation/update, matrix (near) nullspace, and reset functionality */
4071e07b27eSBarry Smith   switch (sr_type) {
4081e07b27eSBarry Smith   case TELESCOPE_DEFAULT:
4091e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
4101e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
4111e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
4121e07b27eSBarry Smith     sred->pctelescope_reset_type              = NULL;
4131e07b27eSBarry Smith     break;
4141e07b27eSBarry Smith   case TELESCOPE_DMDA:
4151e07b27eSBarry Smith     pc->ops->apply                            = PCApply_Telescope_dmda;
416f650675bSDave May     pc->ops->applyrichardson                  = PCApplyRichardson_Telescope_dmda;
4171e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_dmda;
4181e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_dmda;
4191e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda;
4201e07b27eSBarry Smith     sred->pctelescope_reset_type              = PCReset_Telescope_dmda;
4211e07b27eSBarry Smith     break;
422d8b9d5b7SPatrick Sanan   case TELESCOPE_DMPLEX: SETERRQ(comm,PETSC_ERR_SUP,"Support for DMPLEX is currently not available");
4231e07b27eSBarry Smith     break;
424a04a6428SPatrick Sanan   default: SETERRQ(comm,PETSC_ERR_SUP,"Only support for repartitioning DMDA is provided");
4251e07b27eSBarry Smith     break;
4261e07b27eSBarry Smith   }
4271e07b27eSBarry Smith 
4281e07b27eSBarry Smith   /* setup */
4291e07b27eSBarry Smith   if (sred->pctelescope_setup_type) {
4301e07b27eSBarry Smith     ierr = sred->pctelescope_setup_type(pc,sred);CHKERRQ(ierr);
4311e07b27eSBarry Smith   }
4321e07b27eSBarry Smith   /* update */
4331e07b27eSBarry Smith   if (!pc->setupcalled) {
4341e07b27eSBarry Smith     if (sred->pctelescope_matcreate_type) {
4351e07b27eSBarry Smith       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred);CHKERRQ(ierr);
4361e07b27eSBarry Smith     }
4371e07b27eSBarry Smith     if (sred->pctelescope_matnullspacecreate_type) {
438392968a1SPatrick Sanan       ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr);
4391e07b27eSBarry Smith     }
4401e07b27eSBarry Smith   } else {
4411e07b27eSBarry Smith     if (sred->pctelescope_matcreate_type) {
4421e07b27eSBarry Smith       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred);CHKERRQ(ierr);
4431e07b27eSBarry Smith     }
4441e07b27eSBarry Smith   }
4451e07b27eSBarry Smith 
4461e07b27eSBarry Smith   /* common - no construction */
447*2a22aa42SDave May   if (isActiveRank(sred)) {
4481e07b27eSBarry Smith     ierr = KSPSetOperators(sred->ksp,sred->Bred,sred->Bred);CHKERRQ(ierr);
4491e07b27eSBarry Smith     if (pc->setfromoptionscalled && !pc->setupcalled){
4501e07b27eSBarry Smith       ierr = KSPSetFromOptions(sred->ksp);CHKERRQ(ierr);
4511e07b27eSBarry Smith     }
4521e07b27eSBarry Smith   }
4531e07b27eSBarry Smith   PetscFunctionReturn(0);
4541e07b27eSBarry Smith }
4551e07b27eSBarry Smith 
4561e07b27eSBarry Smith static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y)
4571e07b27eSBarry Smith {
4581e07b27eSBarry Smith   PC_Telescope      sred = (PC_Telescope)pc->data;
4591e07b27eSBarry Smith   PetscErrorCode    ierr;
4601e07b27eSBarry Smith   Vec               xtmp,xred,yred;
46113c30530SDave May   PetscInt          i,st,ed;
4621e07b27eSBarry Smith   VecScatter        scatter;
4631e07b27eSBarry Smith   PetscScalar       *array;
4641e07b27eSBarry Smith   const PetscScalar *x_array;
4651e07b27eSBarry Smith 
4661e07b27eSBarry Smith   PetscFunctionBegin;
467bf00f589SPatrick Sanan   ierr = PetscCitationsRegister(citation,&cited);CHKERRQ(ierr);
468bf00f589SPatrick Sanan 
4691e07b27eSBarry Smith   xtmp    = sred->xtmp;
4701e07b27eSBarry Smith   scatter = sred->scatter;
4711e07b27eSBarry Smith   xred    = sred->xred;
4721e07b27eSBarry Smith   yred    = sred->yred;
4731e07b27eSBarry Smith 
4741e07b27eSBarry Smith   /* pull in vector x->xtmp */
4751e07b27eSBarry Smith   ierr = VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4761e07b27eSBarry Smith   ierr = VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4771e07b27eSBarry Smith 
478bf00f589SPatrick Sanan   /* copy vector entries into xred */
4791e07b27eSBarry Smith   ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr);
4801e07b27eSBarry Smith   if (xred) {
4811e07b27eSBarry Smith     PetscScalar *LA_xred;
4821e07b27eSBarry Smith     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
4831e07b27eSBarry Smith     ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr);
4841e07b27eSBarry Smith     for (i=0; i<ed-st; i++) {
4851e07b27eSBarry Smith       LA_xred[i] = x_array[i];
4861e07b27eSBarry Smith     }
4871e07b27eSBarry Smith     ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr);
4881e07b27eSBarry Smith   }
4891e07b27eSBarry Smith   ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr);
4901e07b27eSBarry Smith   /* solve */
491*2a22aa42SDave May   if (isActiveRank(sred)) {
4921e07b27eSBarry Smith     ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr);
493c0decd05SBarry Smith     ierr = KSPCheckSolve(sred->ksp,pc,yred);CHKERRQ(ierr);
4941e07b27eSBarry Smith   }
4951e07b27eSBarry Smith   /* return vector */
4961e07b27eSBarry Smith   ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr);
4971e07b27eSBarry Smith   if (yred) {
4981e07b27eSBarry Smith     const PetscScalar *LA_yred;
4991e07b27eSBarry Smith     ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr);
5001e07b27eSBarry Smith     ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr);
5011e07b27eSBarry Smith     for (i=0; i<ed-st; i++) {
5021e07b27eSBarry Smith       array[i] = LA_yred[i];
5031e07b27eSBarry Smith     }
5041e07b27eSBarry Smith     ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr);
5051e07b27eSBarry Smith   }
5061e07b27eSBarry Smith   ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr);
5071e07b27eSBarry Smith   ierr = VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5081e07b27eSBarry Smith   ierr = VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5091e07b27eSBarry Smith   PetscFunctionReturn(0);
5101e07b27eSBarry Smith }
5111e07b27eSBarry Smith 
512f650675bSDave 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)
513f650675bSDave May {
514f650675bSDave May   PC_Telescope      sred = (PC_Telescope)pc->data;
515f650675bSDave May   PetscErrorCode    ierr;
516a1d91a28SDave May   Vec               xtmp,yred;
517f650675bSDave May   PetscInt          i,st,ed;
518f650675bSDave May   VecScatter        scatter;
519f650675bSDave May   const PetscScalar *x_array;
520f650675bSDave May   PetscBool         default_init_guess_value;
521f650675bSDave May 
522f650675bSDave May   PetscFunctionBegin;
523f650675bSDave May   xtmp    = sred->xtmp;
524f650675bSDave May   scatter = sred->scatter;
525f650675bSDave May   yred    = sred->yred;
526f650675bSDave May 
527f650675bSDave May   if (its > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCApplyRichardson_Telescope only supports max_it = 1");
528f650675bSDave May   *reason = (PCRichardsonConvergedReason)0;
529f650675bSDave May 
530f650675bSDave May   if (!zeroguess) {
531f650675bSDave May     ierr = PetscInfo(pc,"PCTelescope: Scattering y for non-zero initial guess\n");CHKERRQ(ierr);
532f650675bSDave May     /* pull in vector y->xtmp */
533f650675bSDave May     ierr = VecScatterBegin(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
534f650675bSDave May     ierr = VecScatterEnd(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
535f650675bSDave May 
536bf00f589SPatrick Sanan     /* copy vector entries into xred */
537f650675bSDave May     ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr);
538f650675bSDave May     if (yred) {
539f650675bSDave May       PetscScalar *LA_yred;
540f650675bSDave May       ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr);
541f650675bSDave May       ierr = VecGetArray(yred,&LA_yred);CHKERRQ(ierr);
542f650675bSDave May       for (i=0; i<ed-st; i++) {
543f650675bSDave May         LA_yred[i] = x_array[i];
544f650675bSDave May       }
545f650675bSDave May       ierr = VecRestoreArray(yred,&LA_yred);CHKERRQ(ierr);
546f650675bSDave May     }
547f650675bSDave May     ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr);
548f650675bSDave May   }
549f650675bSDave May 
550*2a22aa42SDave May   if (isActiveRank(sred)) {
551f650675bSDave May     ierr = KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value);CHKERRQ(ierr);
552f650675bSDave May     if (!zeroguess) ierr = KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE);CHKERRQ(ierr);
553f650675bSDave May   }
554f650675bSDave May 
555f650675bSDave May   ierr = PCApply_Telescope(pc,x,y);CHKERRQ(ierr);
556f650675bSDave May 
557*2a22aa42SDave May   if (isActiveRank(sred)) {
558f650675bSDave May     ierr = KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value);CHKERRQ(ierr);
559f650675bSDave May   }
560f650675bSDave May 
561f650675bSDave May   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
562f650675bSDave May   *outits = 1;
563f650675bSDave May   PetscFunctionReturn(0);
564f650675bSDave May }
565f650675bSDave May 
5661e07b27eSBarry Smith static PetscErrorCode PCReset_Telescope(PC pc)
5671e07b27eSBarry Smith {
5681e07b27eSBarry Smith   PC_Telescope   sred = (PC_Telescope)pc->data;
5691e07b27eSBarry Smith   PetscErrorCode ierr;
5701e07b27eSBarry Smith 
5711e07b27eSBarry Smith   ierr = ISDestroy(&sred->isin);CHKERRQ(ierr);
5721e07b27eSBarry Smith   ierr = VecScatterDestroy(&sred->scatter);CHKERRQ(ierr);
573e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->xred);CHKERRQ(ierr);
574e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->yred);CHKERRQ(ierr);
575e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->xtmp);CHKERRQ(ierr);
576e3acf2f7SBarry Smith   ierr = MatDestroy(&sred->Bred);CHKERRQ(ierr);
577e3acf2f7SBarry Smith   ierr = KSPReset(sred->ksp);CHKERRQ(ierr);
5781e07b27eSBarry Smith   if (sred->pctelescope_reset_type) {
5791e07b27eSBarry Smith     ierr = sred->pctelescope_reset_type(pc);CHKERRQ(ierr);
5801e07b27eSBarry Smith   }
5811e07b27eSBarry Smith   PetscFunctionReturn(0);
5821e07b27eSBarry Smith }
5831e07b27eSBarry Smith 
5841e07b27eSBarry Smith static PetscErrorCode PCDestroy_Telescope(PC pc)
5851e07b27eSBarry Smith {
5861e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
5871e07b27eSBarry Smith   PetscErrorCode   ierr;
5881e07b27eSBarry Smith 
5891e07b27eSBarry Smith   PetscFunctionBegin;
5901e07b27eSBarry Smith   ierr = PCReset_Telescope(pc);CHKERRQ(ierr);
591e3acf2f7SBarry Smith   ierr = KSPDestroy(&sred->ksp);CHKERRQ(ierr);
5921e07b27eSBarry Smith   ierr = PetscSubcommDestroy(&sred->psubcomm);CHKERRQ(ierr);
593e3acf2f7SBarry Smith   ierr = PetscFree(sred->dm_ctx);CHKERRQ(ierr);
594e3acf2f7SBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
5951e07b27eSBarry Smith   PetscFunctionReturn(0);
5961e07b27eSBarry Smith }
5971e07b27eSBarry Smith 
5984416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Telescope(PetscOptionItems *PetscOptionsObject,PC pc)
5991e07b27eSBarry Smith {
6001e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
6011e07b27eSBarry Smith   PetscErrorCode   ierr;
6021e07b27eSBarry Smith   MPI_Comm         comm;
6031e07b27eSBarry Smith   PetscMPIInt      size;
60448a10b22SPatrick Sanan   PetscBool        flg;
60548a10b22SPatrick Sanan   PetscSubcommType subcommtype;
6061e07b27eSBarry Smith 
6071e07b27eSBarry Smith   PetscFunctionBegin;
6081e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
6091e07b27eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
6101e07b27eSBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Telescope options");CHKERRQ(ierr);
61148a10b22SPatrick Sanan   ierr = PetscOptionsEnum("-pc_telescope_subcomm_type","Subcomm type (interlaced or contiguous)","PCTelescopeSetSubcommType",PetscSubcommTypes,(PetscEnum)sred->subcommtype,(PetscEnum*)&subcommtype,&flg);CHKERRQ(ierr);
61248a10b22SPatrick Sanan   if (flg) {
61348a10b22SPatrick Sanan     ierr = PCTelescopeSetSubcommType(pc,subcommtype);CHKERRQ(ierr);
61448a10b22SPatrick Sanan   }
6151e07b27eSBarry Smith   ierr = PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,0);CHKERRQ(ierr);
6161e07b27eSBarry Smith   if (sred->redfactor > size) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size");
6171e07b27eSBarry Smith   ierr = PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,0);CHKERRQ(ierr);
6187c5279cbSDave May   ierr = PetscOptionsBool("-pc_telescope_ignore_kspcomputeoperators","Ignore method used to compute A","PCTelescopeSetIgnoreKSPComputeOperators",sred->ignore_kspcomputeoperators,&sred->ignore_kspcomputeoperators,0);CHKERRQ(ierr);
6191e07b27eSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
6201e07b27eSBarry Smith   PetscFunctionReturn(0);
6211e07b27eSBarry Smith }
6221e07b27eSBarry Smith 
6231e07b27eSBarry Smith /* PC simplementation specific API's */
6241e07b27eSBarry Smith 
6251e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp)
6261e07b27eSBarry Smith {
6271e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
628bd49479cSSatish Balay   PetscFunctionBegin;
6291e07b27eSBarry Smith   if (ksp) *ksp = red->ksp;
630bd49479cSSatish Balay   PetscFunctionReturn(0);
6311e07b27eSBarry Smith }
6321e07b27eSBarry Smith 
63348a10b22SPatrick Sanan static PetscErrorCode PCTelescopeGetSubcommType_Telescope(PC pc,PetscSubcommType *subcommtype)
63448a10b22SPatrick Sanan {
63548a10b22SPatrick Sanan   PC_Telescope red = (PC_Telescope)pc->data;
63648a10b22SPatrick Sanan   PetscFunctionBegin;
63748a10b22SPatrick Sanan   if (subcommtype) *subcommtype = red->subcommtype;
63848a10b22SPatrick Sanan   PetscFunctionReturn(0);
63948a10b22SPatrick Sanan }
64048a10b22SPatrick Sanan 
64148a10b22SPatrick Sanan static PetscErrorCode PCTelescopeSetSubcommType_Telescope(PC pc,PetscSubcommType subcommtype)
64248a10b22SPatrick Sanan {
64348a10b22SPatrick Sanan   PC_Telescope     red = (PC_Telescope)pc->data;
64448a10b22SPatrick Sanan 
64548a10b22SPatrick Sanan   PetscFunctionBegin;
64648a10b22SPatrick 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.");
64748a10b22SPatrick Sanan   red->subcommtype = subcommtype;
64848a10b22SPatrick Sanan   PetscFunctionReturn(0);
64948a10b22SPatrick Sanan }
65048a10b22SPatrick Sanan 
6511e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact)
6521e07b27eSBarry Smith {
6531e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
654bd49479cSSatish Balay   PetscFunctionBegin;
6551e07b27eSBarry Smith   if (fact) *fact = red->redfactor;
656bd49479cSSatish Balay   PetscFunctionReturn(0);
6571e07b27eSBarry Smith }
6581e07b27eSBarry Smith 
6591e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact)
6601e07b27eSBarry Smith {
6611e07b27eSBarry Smith   PC_Telescope     red = (PC_Telescope)pc->data;
6621e07b27eSBarry Smith   PetscMPIInt      size;
6631e07b27eSBarry Smith   PetscErrorCode   ierr;
6641e07b27eSBarry Smith 
665bd49479cSSatish Balay   PetscFunctionBegin;
6661e07b27eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
6671e07b27eSBarry Smith   if (fact <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact);
6681e07b27eSBarry Smith   if (fact > size) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact);
6691e07b27eSBarry Smith   red->redfactor = fact;
670bd49479cSSatish Balay   PetscFunctionReturn(0);
6711e07b27eSBarry Smith }
6721e07b27eSBarry Smith 
6731e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v)
6741e07b27eSBarry Smith {
6751e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
676bd49479cSSatish Balay   PetscFunctionBegin;
6771e07b27eSBarry Smith   if (v) *v = red->ignore_dm;
678bd49479cSSatish Balay   PetscFunctionReturn(0);
6791e07b27eSBarry Smith }
68048a10b22SPatrick Sanan 
6811e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v)
6821e07b27eSBarry Smith {
6831e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
684bd49479cSSatish Balay   PetscFunctionBegin;
6851e07b27eSBarry Smith   red->ignore_dm = v;
686bd49479cSSatish Balay   PetscFunctionReturn(0);
6871e07b27eSBarry Smith }
6881e07b27eSBarry Smith 
6890ae7c45bSDave May static PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool *v)
6900ae7c45bSDave May {
6910ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
6920ae7c45bSDave May   PetscFunctionBegin;
6930ae7c45bSDave May   if (v) *v = red->ignore_kspcomputeoperators;
6940ae7c45bSDave May   PetscFunctionReturn(0);
6950ae7c45bSDave May }
69648a10b22SPatrick Sanan 
6970ae7c45bSDave May static PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool v)
6980ae7c45bSDave May {
6990ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
7000ae7c45bSDave May   PetscFunctionBegin;
7010ae7c45bSDave May   red->ignore_kspcomputeoperators = v;
7020ae7c45bSDave May   PetscFunctionReturn(0);
7030ae7c45bSDave May }
7040ae7c45bSDave May 
7051e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm)
7061e07b27eSBarry Smith {
7071e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
708bd49479cSSatish Balay   PetscFunctionBegin;
7091e07b27eSBarry Smith   *dm = private_PCTelescopeGetSubDM(red);
710bd49479cSSatish Balay   PetscFunctionReturn(0);
7111e07b27eSBarry Smith }
7121e07b27eSBarry Smith 
7131e07b27eSBarry Smith /*@
7141e07b27eSBarry Smith  PCTelescopeGetKSP - Gets the KSP created by the telescoping PC.
7151e07b27eSBarry Smith 
7161e07b27eSBarry Smith  Not Collective
7171e07b27eSBarry Smith 
7181e07b27eSBarry Smith  Input Parameter:
7191e07b27eSBarry Smith .  pc - the preconditioner context
7201e07b27eSBarry Smith 
7211e07b27eSBarry Smith  Output Parameter:
7221e07b27eSBarry Smith .  subksp - the KSP defined the smaller set of processes
7231e07b27eSBarry Smith 
7241e07b27eSBarry Smith  Level: advanced
7251e07b27eSBarry Smith 
7261e07b27eSBarry Smith .keywords: PC, telescopting solve
7271e07b27eSBarry Smith @*/
7281e07b27eSBarry Smith PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp)
7291e07b27eSBarry Smith {
730bd49479cSSatish Balay   PetscErrorCode ierr;
731bd49479cSSatish Balay   PetscFunctionBegin;
732163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp));CHKERRQ(ierr);
733bd49479cSSatish Balay   PetscFunctionReturn(0);
7341e07b27eSBarry Smith }
7351e07b27eSBarry Smith 
7361e07b27eSBarry Smith /*@
7371e07b27eSBarry Smith  PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by.
7381e07b27eSBarry Smith 
7391e07b27eSBarry Smith  Not Collective
7401e07b27eSBarry Smith 
7411e07b27eSBarry Smith  Input Parameter:
7421e07b27eSBarry Smith .  pc - the preconditioner context
7431e07b27eSBarry Smith 
7441e07b27eSBarry Smith  Output Parameter:
7451e07b27eSBarry Smith .  fact - the reduction factor
7461e07b27eSBarry Smith 
7471e07b27eSBarry Smith  Level: advanced
7481e07b27eSBarry Smith 
7491e07b27eSBarry Smith .keywords: PC, telescoping solve
7501e07b27eSBarry Smith @*/
7511e07b27eSBarry Smith PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact)
7521e07b27eSBarry Smith {
753bd49479cSSatish Balay   PetscErrorCode ierr;
754bd49479cSSatish Balay   PetscFunctionBegin;
755163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact));CHKERRQ(ierr);
756bd49479cSSatish Balay   PetscFunctionReturn(0);
7571e07b27eSBarry Smith }
7581e07b27eSBarry Smith 
7591e07b27eSBarry Smith /*@
7601e07b27eSBarry Smith  PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by.
7611e07b27eSBarry Smith 
7621e07b27eSBarry Smith  Not Collective
7631e07b27eSBarry Smith 
7641e07b27eSBarry Smith  Input Parameter:
7651e07b27eSBarry Smith .  pc - the preconditioner context
7661e07b27eSBarry Smith 
7671e07b27eSBarry Smith  Output Parameter:
7681e07b27eSBarry Smith .  fact - the reduction factor
7691e07b27eSBarry Smith 
7701e07b27eSBarry Smith  Level: advanced
7711e07b27eSBarry Smith 
7721e07b27eSBarry Smith .keywords: PC, telescoping solve
7731e07b27eSBarry Smith @*/
7741e07b27eSBarry Smith PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact)
7751e07b27eSBarry Smith {
776bd49479cSSatish Balay   PetscErrorCode ierr;
777bd49479cSSatish Balay   PetscFunctionBegin;
778bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact));CHKERRQ(ierr);
779bd49479cSSatish Balay   PetscFunctionReturn(0);
7801e07b27eSBarry Smith }
7811e07b27eSBarry Smith 
7821e07b27eSBarry Smith /*@
7831e07b27eSBarry Smith  PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used.
7841e07b27eSBarry Smith 
7851e07b27eSBarry Smith  Not Collective
7861e07b27eSBarry Smith 
7871e07b27eSBarry Smith  Input Parameter:
7881e07b27eSBarry Smith .  pc - the preconditioner context
7891e07b27eSBarry Smith 
7901e07b27eSBarry Smith  Output Parameter:
7911e07b27eSBarry Smith .  v - the flag
7921e07b27eSBarry Smith 
7931e07b27eSBarry Smith  Level: advanced
7941e07b27eSBarry Smith 
7951e07b27eSBarry Smith .keywords: PC, telescoping solve
7961e07b27eSBarry Smith @*/
7971e07b27eSBarry Smith PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v)
7981e07b27eSBarry Smith {
799bd49479cSSatish Balay   PetscErrorCode ierr;
800bd49479cSSatish Balay   PetscFunctionBegin;
801163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr);
802bd49479cSSatish Balay   PetscFunctionReturn(0);
8031e07b27eSBarry Smith }
8041e07b27eSBarry Smith 
8051e07b27eSBarry Smith /*@
8061e07b27eSBarry Smith  PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC.
8071e07b27eSBarry Smith 
8081e07b27eSBarry Smith  Not Collective
8091e07b27eSBarry Smith 
8101e07b27eSBarry Smith  Input Parameter:
8111e07b27eSBarry Smith .  pc - the preconditioner context
8121e07b27eSBarry Smith 
8131e07b27eSBarry Smith  Output Parameter:
8141e07b27eSBarry Smith .  v - Use PETSC_TRUE to ignore any DM
8151e07b27eSBarry Smith 
8161e07b27eSBarry Smith  Level: advanced
8171e07b27eSBarry Smith 
8181e07b27eSBarry Smith .keywords: PC, telescoping solve
8191e07b27eSBarry Smith @*/
820bfd6bcc6SSatish Balay PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v)
8211e07b27eSBarry Smith {
822bd49479cSSatish Balay   PetscErrorCode ierr;
823bd49479cSSatish Balay   PetscFunctionBegin;
824bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr);
825bd49479cSSatish Balay   PetscFunctionReturn(0);
8261e07b27eSBarry Smith }
8271e07b27eSBarry Smith 
8281e07b27eSBarry Smith /*@
8290ae7c45bSDave May  PCTelescopeGetIgnoreKSPComputeOperators - Get the flag indicating if KSPComputeOperators will be used.
8300ae7c45bSDave May 
8310ae7c45bSDave May  Not Collective
8320ae7c45bSDave May 
8330ae7c45bSDave May  Input Parameter:
8340ae7c45bSDave May .  pc - the preconditioner context
8350ae7c45bSDave May 
8360ae7c45bSDave May  Output Parameter:
8370ae7c45bSDave May .  v - the flag
8380ae7c45bSDave May 
8390ae7c45bSDave May  Level: advanced
8400ae7c45bSDave May 
8410ae7c45bSDave May .keywords: PC, telescoping solve
8420ae7c45bSDave May @*/
8430ae7c45bSDave May PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC pc,PetscBool *v)
8440ae7c45bSDave May {
8450ae7c45bSDave May   PetscErrorCode ierr;
8460ae7c45bSDave May   PetscFunctionBegin;
847163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr);
8480ae7c45bSDave May   PetscFunctionReturn(0);
8490ae7c45bSDave May }
8500ae7c45bSDave May 
8510ae7c45bSDave May /*@
8520ae7c45bSDave May  PCTelescopeSetIgnoreKSPComputeOperators - Set a flag to ignore KSPComputeOperators.
8530ae7c45bSDave May 
8540ae7c45bSDave May  Not Collective
8550ae7c45bSDave May 
8560ae7c45bSDave May  Input Parameter:
8570ae7c45bSDave May .  pc - the preconditioner context
8580ae7c45bSDave May 
8590ae7c45bSDave May  Output Parameter:
860a954d8f4SDave May .  v - Use PETSC_TRUE to ignore the method (if defined) set via KSPSetComputeOperators on pc
8610ae7c45bSDave May 
8620ae7c45bSDave May  Level: advanced
8630ae7c45bSDave May 
8640ae7c45bSDave May .keywords: PC, telescoping solve
8650ae7c45bSDave May @*/
8660ae7c45bSDave May PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC pc,PetscBool v)
8670ae7c45bSDave May {
8680ae7c45bSDave May   PetscErrorCode ierr;
8690ae7c45bSDave May   PetscFunctionBegin;
8700ae7c45bSDave May   ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr);
8710ae7c45bSDave May   PetscFunctionReturn(0);
8720ae7c45bSDave May }
8730ae7c45bSDave May 
8740ae7c45bSDave May /*@
8751e07b27eSBarry Smith  PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP.
8761e07b27eSBarry Smith 
8771e07b27eSBarry Smith  Not Collective
8781e07b27eSBarry Smith 
8791e07b27eSBarry Smith  Input Parameter:
8801e07b27eSBarry Smith .  pc - the preconditioner context
8811e07b27eSBarry Smith 
8821e07b27eSBarry Smith  Output Parameter:
8831e07b27eSBarry Smith .  subdm - The re-partitioned DM
8841e07b27eSBarry Smith 
8851e07b27eSBarry Smith  Level: advanced
8861e07b27eSBarry Smith 
8871e07b27eSBarry Smith .keywords: PC, telescoping solve
8881e07b27eSBarry Smith @*/
8891e07b27eSBarry Smith PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm)
8901e07b27eSBarry Smith {
891bd49479cSSatish Balay   PetscErrorCode ierr;
892bd49479cSSatish Balay   PetscFunctionBegin;
893163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm));CHKERRQ(ierr);
894bd49479cSSatish Balay   PetscFunctionReturn(0);
8951e07b27eSBarry Smith }
8961e07b27eSBarry Smith 
89748a10b22SPatrick Sanan /*@
89848a10b22SPatrick Sanan  PCTelescopeSetSubcommType - set subcommunicator type (interlaced or contiguous)
89948a10b22SPatrick Sanan 
90048a10b22SPatrick Sanan  Logically Collective
90148a10b22SPatrick Sanan 
90248a10b22SPatrick Sanan  Input Parameter:
9031dae98e4SBarry Smith +  pc - the preconditioner context
9041dae98e4SBarry Smith -  subcommtype - the subcommunicator type (see PetscSubcommType)
90548a10b22SPatrick Sanan 
90648a10b22SPatrick Sanan  Level: advanced
90748a10b22SPatrick Sanan 
90848a10b22SPatrick Sanan .keywords: PC, telescoping solve
90948a10b22SPatrick Sanan 
91048a10b22SPatrick Sanan .seealso: PetscSubcommType, PetscSubcomm, PCTELESCOPE
91148a10b22SPatrick Sanan @*/
91248a10b22SPatrick Sanan PetscErrorCode PCTelescopeSetSubcommType(PC pc, PetscSubcommType subcommtype)
91348a10b22SPatrick Sanan {
91448a10b22SPatrick Sanan   PetscErrorCode ierr;
91548a10b22SPatrick Sanan   PetscFunctionBegin;
91648a10b22SPatrick Sanan   ierr = PetscTryMethod(pc,"PCTelescopeSetSubcommType_C",(PC,PetscSubcommType),(pc,subcommtype));CHKERRQ(ierr);
91748a10b22SPatrick Sanan   PetscFunctionReturn(0);
91848a10b22SPatrick Sanan }
91948a10b22SPatrick Sanan 
92048a10b22SPatrick Sanan /*@
92148a10b22SPatrick Sanan  PCTelescopeGetSubcommType - Get the subcommunicator type (interlaced or contiguous)
92248a10b22SPatrick Sanan 
92348a10b22SPatrick Sanan  Not Collective
92448a10b22SPatrick Sanan 
92548a10b22SPatrick Sanan  Input Parameter:
92648a10b22SPatrick Sanan .  pc - the preconditioner context
92748a10b22SPatrick Sanan 
92848a10b22SPatrick Sanan  Output Parameter:
92948a10b22SPatrick Sanan .  subcommtype - the subcommunicator type (see PetscSubcommType)
93048a10b22SPatrick Sanan 
93148a10b22SPatrick Sanan  Level: advanced
93248a10b22SPatrick Sanan 
93348a10b22SPatrick Sanan .keywords: PC, telescoping solve
93448a10b22SPatrick Sanan 
9351dae98e4SBarry Smith .seealso: PetscSubcomm, PetscSubcommType, PCTELESCOPE
93648a10b22SPatrick Sanan @*/
9371dae98e4SBarry Smith PetscErrorCode PCTelescopeGetSubcommType(PC pc, PetscSubcommType *subcommtype)
93848a10b22SPatrick Sanan {
93948a10b22SPatrick Sanan   PetscErrorCode ierr;
94048a10b22SPatrick Sanan   PetscFunctionBegin;
94148a10b22SPatrick Sanan   ierr = PetscUseMethod(pc,"PCTelescopeGetSubcommType_C",(PC,PetscSubcommType*),(pc,subcommtype));CHKERRQ(ierr);
94248a10b22SPatrick Sanan   PetscFunctionReturn(0);
94348a10b22SPatrick Sanan }
94448a10b22SPatrick Sanan 
9451e07b27eSBarry Smith /* -------------------------------------------------------------------------------------*/
9461e07b27eSBarry Smith /*MC
9471e07b27eSBarry Smith    PCTELESCOPE - Runs a KSP solver on a sub-group of processors. MPI processes not in the sub-communicator are idle during the solve.
9481e07b27eSBarry Smith 
9491e07b27eSBarry Smith    Options Database:
950a04a6428SPatrick Sanan +  -pc_telescope_reduction_factor <r> - factor to use communicator size by. e.g. with 64 MPI processes and r=4, the new sub-communicator will have 64/4 = 16 ranks.
951a04a6428SPatrick Sanan -  -pc_telescope_ignore_dm  - flag to indicate whether an attached DM should be ignored
952a04a6428SPatrick Sanan -  -pc_telescope_subcomm_type <interlaced,contiguous> - how to define the reduced communicator. see PetscSubcomm for more.
9531e07b27eSBarry Smith 
9541e07b27eSBarry Smith    Level: advanced
9551e07b27eSBarry Smith 
9561e07b27eSBarry Smith    Notes:
9576fc41876SBarry Smith    The preconditioner is deemed telescopic as it only calls KSPSolve() on a single
9588439623fSPatrick Sanan    sub-communicator, in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators.
9598439623fSPatrick Sanan    This means there will be MPI processes which will be idle during the application of this preconditioner.
9606fc41876SBarry Smith 
9611e07b27eSBarry Smith    The default KSP is PREONLY. If a DM is attached to the PC, it is re-partitioned on the sub-communicator.
9621e07b27eSBarry Smith    Both the Bmat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned DM.
9631e07b27eSBarry Smith    Currently only support for re-partitioning a DMDA is provided.
9648439623fSPatrick Sanan    Any nullspace attached to the original Bmat operator is extracted, re-partitioned and set on the repartitioned Bmat operator.
9651e07b27eSBarry Smith    KSPSetComputeOperators() is not propagated to the sub KSP.
9661e07b27eSBarry Smith    Currently there is no support for the flag -pc_use_amat
9671e07b27eSBarry Smith 
9686fc41876SBarry Smith    Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation
9698439623fSPatrick Sanan    creates a child sub-communicator (c') containing fewer MPI processes than the original parent preconditioner (PC).
9706fc41876SBarry Smith 
9716fc41876SBarry Smith   Developer Notes:
9726fc41876SBarry Smith    During PCSetup, the B operator is scattered onto c'.
9736fc41876SBarry Smith    Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c').
9748439623fSPatrick Sanan    Then, KSPSolve() is executed on the c' communicator.
9756fc41876SBarry Smith 
9766fc41876SBarry Smith    The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED
977a04a6428SPatrick 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.
9786fc41876SBarry Smith 
979005d9f20SPatrick Sanan    The telescoping preconditioner is aware of nullspaces and near nullspaces which are attached to the B operator.
9808439623fSPatrick Sanan    In the case where B has a (near) nullspace attached, the (near) nullspace vectors are extracted from B and mapped into
981005d9f20SPatrick Sanan    a new (near) nullspace, defined on the sub-communicator, which is attached to B' (the B operator which was scattered to c')
9826fc41876SBarry Smith 
9836fc41876SBarry Smith    The telescoping preconditioner is aware of an attached DM. In the event that the DM is of type DMDA (2D or 3D -
9846fc41876SBarry Smith    1D support for 1D DMDAs is not provided), a new DMDA is created on c' (e.g. it is re-partitioned), and this new DM
9856fc41876SBarry Smith    is attached the sub KSPSolve(). The design of telescope is such that it should be possible to extend support
9868439623fSPatrick Sanan    for re-partitioning other to DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs.
9876fc41876SBarry Smith 
9886fc41876SBarry Smith    By default, B' is defined by simply fusing rows from different MPI processes
9896fc41876SBarry Smith 
9908439623fSPatrick Sanan    When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permutation of B
9917dae84e0SHong Zhang    into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatCreateSubMatrices(), (iii) fusing the
9926fc41876SBarry Smith    locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat()
9936fc41876SBarry Smith 
9948439623fSPatrick Sanan    Limitations/improvements include the following.
9958439623fSPatrick Sanan    VecPlaceArray() could be used within PCApply() to improve efficiency and reduce memory usage.
9966fc41876SBarry Smith 
9976fc41876SBarry Smith    The symmetric permutation used when a DMDA is encountered is performed via explicitly assmbleming a permutation matrix P,
9988439623fSPatrick 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
9998439623fSPatrick Sanan    VecPermute() does not supported for the use case required here. By computing P, one can permute both the operator and RHS in a
10006fc41876SBarry Smith    consistent manner.
10016fc41876SBarry Smith 
10028439623fSPatrick Sanan    Mapping of vectors is performed in the following way.
10038439623fSPatrick Sanan    Suppose the parent comm size was 4, and we set a reduction factor of 2; this would give a comm size on c' of 2.
10048439623fSPatrick Sanan    Using the interlaced creation routine, the ranks in c with color = 0 will be rank 0 and 2.
10058439623fSPatrick Sanan    We perform the scatter to the sub-comm in the following way.
10066fc41876SBarry Smith    [1] Given a vector x defined on comm c
10076fc41876SBarry Smith 
10086fc41876SBarry Smith    rank(c) : _________ 0 ______  ________ 1 _______  ________ 2 _____________ ___________ 3 __________
10096fc41876SBarry Smith          x : [0, 1, 2, 3, 4, 5] [6, 7, 8, 9, 10, 11] [12, 13, 14, 15, 16, 17] [18, 19, 20, 21, 22, 23]
10106fc41876SBarry Smith 
10116fc41876SBarry Smith    scatter to xtmp defined also on comm c so that we have the following values
10126fc41876SBarry Smith 
10136fc41876SBarry Smith    rank(c) : ___________________ 0 ________________  _1_  ______________________ 2 _______________________  __3_
10146fc41876SBarry Smith       xtmp : [0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11] [  ] [12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23] [  ]
10156fc41876SBarry Smith 
10166fc41876SBarry Smith    The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') have no values
10176fc41876SBarry Smith 
10186fc41876SBarry Smith 
10196fc41876SBarry Smith    [2] Copy the value from rank 0, 2 (indices with respect to comm c) into the vector xred which is defined on communicator c'.
10206fc41876SBarry Smith    Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0.
10216fc41876SBarry Smith 
10226fc41876SBarry Smith     rank(c') : ___________________ 0 _______________  ______________________ 1 _____________________
10236fc41876SBarry Smith       xred : [0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11] [12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23]
10246fc41876SBarry Smith 
10251e07b27eSBarry Smith 
10261e07b27eSBarry Smith   Contributed by Dave May
10271e07b27eSBarry Smith 
1028bf00f589SPatrick Sanan   Reference:
1029bf00f589SPatrick 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
1030bf00f589SPatrick Sanan 
10316fc41876SBarry Smith .seealso:  PCTelescopeGetKSP(), PCTelescopeGetDM(), PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(), PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT
10321e07b27eSBarry Smith M*/
10331e07b27eSBarry Smith PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc)
10341e07b27eSBarry Smith {
10351e07b27eSBarry Smith   PetscErrorCode       ierr;
10361e07b27eSBarry Smith   struct _PC_Telescope *sred;
10371e07b27eSBarry Smith 
10381e07b27eSBarry Smith   PetscFunctionBegin;
10391e07b27eSBarry Smith   ierr = PetscNewLog(pc,&sred);CHKERRQ(ierr);
1040*2a22aa42SDave May   sred->psubcomm       = NULL;
104148a10b22SPatrick Sanan   sred->subcommtype    = PETSC_SUBCOMM_INTERLACED;
1042*2a22aa42SDave May   sred->subcomm        = MPI_COMM_NULL;
10431e07b27eSBarry Smith   sred->redfactor      = 1;
10441e07b27eSBarry Smith   sred->ignore_dm      = PETSC_FALSE;
10457c5279cbSDave May   sred->ignore_kspcomputeoperators = PETSC_FALSE;
10461e07b27eSBarry Smith   pc->data             = (void*)sred;
10471e07b27eSBarry Smith 
10481e07b27eSBarry Smith   pc->ops->apply           = PCApply_Telescope;
10491e07b27eSBarry Smith   pc->ops->applytranspose  = NULL;
1050f650675bSDave May   pc->ops->applyrichardson = PCApplyRichardson_Telescope;
10511e07b27eSBarry Smith   pc->ops->setup           = PCSetUp_Telescope;
10521e07b27eSBarry Smith   pc->ops->destroy         = PCDestroy_Telescope;
10531e07b27eSBarry Smith   pc->ops->reset           = PCReset_Telescope;
10541e07b27eSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Telescope;
10551e07b27eSBarry Smith   pc->ops->view            = PCView_Telescope;
10561e07b27eSBarry Smith 
10571e07b27eSBarry Smith   sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
10581e07b27eSBarry Smith   sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
10591e07b27eSBarry Smith   sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
10601e07b27eSBarry Smith   sred->pctelescope_reset_type              = NULL;
10611e07b27eSBarry Smith 
10621e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetKSP_C",PCTelescopeGetKSP_Telescope);CHKERRQ(ierr);
106348a10b22SPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetSubcommType_C",PCTelescopeGetSubcommType_Telescope);CHKERRQ(ierr);
106448a10b22SPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetSubcommType_C",PCTelescopeSetSubcommType_Telescope);CHKERRQ(ierr);
10651e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetReductionFactor_C",PCTelescopeGetReductionFactor_Telescope);CHKERRQ(ierr);
10661e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetReductionFactor_C",PCTelescopeSetReductionFactor_Telescope);CHKERRQ(ierr);
10671e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreDM_C",PCTelescopeGetIgnoreDM_Telescope);CHKERRQ(ierr);
10681e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreDM_C",PCTelescopeSetIgnoreDM_Telescope);CHKERRQ(ierr);
10690ae7c45bSDave May   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",PCTelescopeGetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr);
10700ae7c45bSDave May   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",PCTelescopeSetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr);
10711e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetDM_C",PCTelescopeGetDM_Telescope);CHKERRQ(ierr);
10721e07b27eSBarry Smith   PetscFunctionReturn(0);
10731e07b27eSBarry Smith }
1074