xref: /petsc/src/ksp/pc/impls/telescope/telescope.c (revision 48a10b22724dcb2133e337caacd0d29403d7be99)
11e07b27eSBarry Smith 
21e07b27eSBarry Smith 
31e07b27eSBarry Smith 
4120bdd93SDave May #include <petsc/private/matimpl.h>
56fc41876SBarry Smith #include <petsc/private/pcimpl.h>
61e07b27eSBarry Smith #include <petscksp.h> /*I "petscksp.h" I*/
71e07b27eSBarry Smith #include <petscdm.h> /*I "petscdm.h" I*/
8*48a10b22SPatrick Sanan #include <petscsys.h>
91e07b27eSBarry Smith 
101e07b27eSBarry Smith #include "telescope.h"
111e07b27eSBarry Smith 
121e07b27eSBarry Smith /*
131e07b27eSBarry Smith  PCTelescopeSetUp_default()
141e07b27eSBarry Smith  PCTelescopeMatCreate_default()
151e07b27eSBarry Smith 
161e07b27eSBarry Smith  default
171e07b27eSBarry Smith 
181e07b27eSBarry Smith  // scatter in
191e07b27eSBarry Smith  x(comm) -> xtmp(comm)
201e07b27eSBarry Smith 
211e07b27eSBarry Smith  xred(subcomm) <- xtmp
221e07b27eSBarry Smith  yred(subcomm)
231e07b27eSBarry Smith 
241e07b27eSBarry Smith  yred(subcomm) --> xtmp
251e07b27eSBarry Smith 
261e07b27eSBarry Smith  // scatter out
271e07b27eSBarry Smith  xtmp(comm) -> y(comm)
281e07b27eSBarry Smith */
291e07b27eSBarry Smith 
301e07b27eSBarry Smith PetscBool isActiveRank(PetscSubcomm scomm)
311e07b27eSBarry Smith {
321e07b27eSBarry Smith   if (scomm->color == 0) { return PETSC_TRUE; }
331e07b27eSBarry Smith   else { return PETSC_FALSE; }
341e07b27eSBarry Smith }
351e07b27eSBarry Smith 
361e07b27eSBarry Smith #undef __FUNCT__
371e07b27eSBarry Smith #define __FUNCT__ "private_PCTelescopeGetSubDM"
381e07b27eSBarry Smith DM private_PCTelescopeGetSubDM(PC_Telescope sred)
391e07b27eSBarry Smith {
40c6a0d831SBarry Smith   DM subdm = NULL;
411e07b27eSBarry Smith 
421e07b27eSBarry Smith   if (!isActiveRank(sred->psubcomm)) { subdm = NULL; }
431e07b27eSBarry Smith   else {
441e07b27eSBarry Smith     switch (sred->sr_type) {
451e07b27eSBarry Smith     case TELESCOPE_DEFAULT: subdm = NULL;
461e07b27eSBarry Smith       break;
471e07b27eSBarry Smith     case TELESCOPE_DMDA:    subdm = ((PC_Telescope_DMDACtx*)sred->dm_ctx)->dmrepart;
481e07b27eSBarry Smith       break;
491e07b27eSBarry Smith     case TELESCOPE_DMPLEX:  subdm = NULL;
501e07b27eSBarry Smith       break;
511e07b27eSBarry Smith     }
521e07b27eSBarry Smith   }
531e07b27eSBarry Smith   return(subdm);
541e07b27eSBarry Smith }
551e07b27eSBarry Smith 
561e07b27eSBarry Smith #undef __FUNCT__
571e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeSetUp_default"
581e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_default(PC pc,PC_Telescope sred)
591e07b27eSBarry Smith {
601e07b27eSBarry Smith   PetscErrorCode ierr;
611e07b27eSBarry Smith   PetscInt       m,M,bs,st,ed;
621e07b27eSBarry Smith   Vec            x,xred,yred,xtmp;
631e07b27eSBarry Smith   Mat            B;
641e07b27eSBarry Smith   MPI_Comm       comm,subcomm;
651e07b27eSBarry Smith   VecScatter     scatter;
661e07b27eSBarry Smith   IS             isin;
671e07b27eSBarry Smith 
681e07b27eSBarry Smith   PetscFunctionBegin;
691e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: setup (default)\n");CHKERRQ(ierr);
701e07b27eSBarry Smith   comm = PetscSubcommParent(sred->psubcomm);
711e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
721e07b27eSBarry Smith 
731e07b27eSBarry Smith   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
741e07b27eSBarry Smith   ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr);
751e07b27eSBarry Smith   ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr);
761e07b27eSBarry Smith   ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr);
771e07b27eSBarry Smith 
781e07b27eSBarry Smith   xred = NULL;
793ac26c5eSBarry Smith   m    = 0;
801e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
811e07b27eSBarry Smith     ierr = VecCreate(subcomm,&xred);CHKERRQ(ierr);
821e07b27eSBarry Smith     ierr = VecSetSizes(xred,PETSC_DECIDE,M);CHKERRQ(ierr);
831e07b27eSBarry Smith     ierr = VecSetBlockSize(xred,bs);CHKERRQ(ierr);
841e07b27eSBarry Smith     ierr = VecSetFromOptions(xred);CHKERRQ(ierr);
85ca43db0aSBarry Smith     ierr = VecGetLocalSize(xred,&m);CHKERRQ(ierr);
861e07b27eSBarry Smith   }
871e07b27eSBarry Smith 
881e07b27eSBarry Smith   yred = NULL;
891e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
901e07b27eSBarry Smith     ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr);
911e07b27eSBarry Smith   }
921e07b27eSBarry Smith 
931e07b27eSBarry Smith   ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr);
941e07b27eSBarry Smith   ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr);
951e07b27eSBarry Smith   ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr);
961e07b27eSBarry Smith   ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr);
971e07b27eSBarry Smith 
981e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
991e07b27eSBarry Smith     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
1001e07b27eSBarry Smith     ierr = ISCreateStride(comm,(ed-st),st,1,&isin);CHKERRQ(ierr);
1011e07b27eSBarry Smith   } else {
1021e07b27eSBarry Smith     ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr);
1033ac26c5eSBarry Smith     ierr = ISCreateStride(comm,0,st,1,&isin);CHKERRQ(ierr);
1041e07b27eSBarry Smith   }
1051e07b27eSBarry Smith   ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr);
1061e07b27eSBarry Smith 
1071e07b27eSBarry Smith   ierr = VecScatterCreate(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr);
1081e07b27eSBarry Smith 
1091e07b27eSBarry Smith   sred->isin    = isin;
1101e07b27eSBarry Smith   sred->scatter = scatter;
1111e07b27eSBarry Smith   sred->xred    = xred;
1121e07b27eSBarry Smith   sred->yred    = yred;
1131e07b27eSBarry Smith   sred->xtmp    = xtmp;
1141e07b27eSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
1151e07b27eSBarry Smith   PetscFunctionReturn(0);
1161e07b27eSBarry Smith }
1171e07b27eSBarry Smith 
1181e07b27eSBarry Smith #undef __FUNCT__
1191e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeMatCreate_default"
1201e07b27eSBarry Smith PetscErrorCode PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A)
1211e07b27eSBarry Smith {
1221e07b27eSBarry Smith   PetscErrorCode ierr;
1231e07b27eSBarry Smith   MPI_Comm       comm,subcomm;
1241e07b27eSBarry Smith   Mat            Bred,B;
1251e07b27eSBarry Smith   PetscInt       nr,nc;
1261e07b27eSBarry Smith   IS             isrow,iscol;
1271e07b27eSBarry Smith   Mat            Blocal,*_Blocal;
1281e07b27eSBarry Smith 
1291e07b27eSBarry Smith   PetscFunctionBegin;
1301e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (default)\n");CHKERRQ(ierr);
1311e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1321e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
1331e07b27eSBarry Smith   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
1341e07b27eSBarry Smith   ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr);
1351e07b27eSBarry Smith   isrow = sred->isin;
1361e07b27eSBarry Smith   ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr);
1371e07b27eSBarry Smith   ierr = MatGetSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr);
1381e07b27eSBarry Smith   Blocal = *_Blocal;
1391e07b27eSBarry Smith   ierr = PetscFree(_Blocal);CHKERRQ(ierr);
1401e07b27eSBarry Smith   Bred = NULL;
1411e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
1421e07b27eSBarry Smith     PetscInt mm;
1431e07b27eSBarry Smith 
1441e07b27eSBarry Smith     if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; }
1451e07b27eSBarry Smith 
1461e07b27eSBarry Smith     ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr);
1471e07b27eSBarry Smith     ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr);
1481e07b27eSBarry Smith   }
1491e07b27eSBarry Smith   *A = Bred;
1501e07b27eSBarry Smith   ierr = ISDestroy(&iscol);CHKERRQ(ierr);
1511e07b27eSBarry Smith   ierr = MatDestroy(&Blocal);CHKERRQ(ierr);
1521e07b27eSBarry Smith   PetscFunctionReturn(0);
1531e07b27eSBarry Smith }
1541e07b27eSBarry Smith 
1551e07b27eSBarry Smith #undef __FUNCT__
1561e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeMatNullSpaceCreate_default"
1571e07b27eSBarry Smith PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat)
1581e07b27eSBarry Smith {
1591e07b27eSBarry Smith   PetscErrorCode   ierr;
1601e07b27eSBarry Smith   MatNullSpace     nullspace,sub_nullspace;
1611e07b27eSBarry Smith   Mat              A,B;
1621e07b27eSBarry Smith   PetscBool        has_const;
163a947c41eSDave May   PetscInt         i,k,n = 0;
1641e07b27eSBarry Smith   const Vec        *vecs;
165c41e779fSDave May   Vec              *sub_vecs = NULL;
1661e07b27eSBarry Smith   MPI_Comm         subcomm;
1671e07b27eSBarry Smith 
1681e07b27eSBarry Smith   PetscFunctionBegin;
1691e07b27eSBarry Smith   ierr = PCGetOperators(pc,&A,&B);CHKERRQ(ierr);
1701e07b27eSBarry Smith   ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr);
171a947c41eSDave May   if (!nullspace) PetscFunctionReturn(0);
1721e07b27eSBarry Smith 
1731e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: generating nullspace (default)\n");CHKERRQ(ierr);
1741e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
1751e07b27eSBarry Smith   ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr);
1761e07b27eSBarry Smith 
1771e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
178e3acf2f7SBarry Smith     if (n) {
179e3acf2f7SBarry Smith       ierr = VecDuplicateVecs(sred->xred,n,&sub_vecs);CHKERRQ(ierr);
1801e07b27eSBarry Smith     }
1811e07b27eSBarry Smith   }
1821e07b27eSBarry Smith 
1831e07b27eSBarry Smith   /* copy entries */
1841e07b27eSBarry Smith   for (k=0; k<n; k++) {
1851e07b27eSBarry Smith     const PetscScalar *x_array;
1861e07b27eSBarry Smith     PetscScalar       *LA_sub_vec;
18713c30530SDave May     PetscInt          st,ed;
1881e07b27eSBarry Smith 
1891e07b27eSBarry Smith     /* pull in vector x->xtmp */
1901e07b27eSBarry Smith     ierr = VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1911e07b27eSBarry Smith     ierr = VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19247856c66SBarry Smith     if (sub_vecs) {
1931e07b27eSBarry Smith       /* copy vector entires into xred */
1941e07b27eSBarry Smith       ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
195ea2b237eSDave May       if (sub_vecs[k]) {
1961e07b27eSBarry Smith         ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr);
1971e07b27eSBarry Smith         ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
1981e07b27eSBarry Smith         for (i=0; i<ed-st; i++) {
1991e07b27eSBarry Smith           LA_sub_vec[i] = x_array[i];
2001e07b27eSBarry Smith         }
2011e07b27eSBarry Smith         ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
2021e07b27eSBarry Smith       }
2031e07b27eSBarry Smith       ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
2041e07b27eSBarry Smith     }
20547856c66SBarry Smith   }
2061e07b27eSBarry Smith 
2071e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
2081e07b27eSBarry Smith     /* create new nullspace for redundant object */
2091e07b27eSBarry Smith     ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,&sub_nullspace);CHKERRQ(ierr);
210120bdd93SDave May     sub_nullspace->remove = nullspace->remove;
211120bdd93SDave May     sub_nullspace->rmctx = nullspace->rmctx;
212120bdd93SDave May 
2131e07b27eSBarry Smith     /* attach redundant nullspace to Bred */
2141e07b27eSBarry Smith     ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr);
215e3acf2f7SBarry Smith     ierr = VecDestroyVecs(n,&sub_vecs);CHKERRQ(ierr);
2161e07b27eSBarry Smith   }
2171e07b27eSBarry Smith   PetscFunctionReturn(0);
2181e07b27eSBarry Smith }
2191e07b27eSBarry Smith 
2201e07b27eSBarry Smith #undef __FUNCT__
2211e07b27eSBarry Smith #define __FUNCT__ "PCView_Telescope"
2221e07b27eSBarry Smith static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer)
2231e07b27eSBarry Smith {
2241e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
2251e07b27eSBarry Smith   PetscErrorCode   ierr;
2261e07b27eSBarry Smith   PetscBool        iascii,isstring;
2271e07b27eSBarry Smith   PetscViewer      subviewer;
2281e07b27eSBarry Smith 
2291e07b27eSBarry Smith   PetscFunctionBegin;
2301e07b27eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
2311e07b27eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
2321e07b27eSBarry Smith   if (iascii) {
2331e07b27eSBarry Smith     if (!sred->psubcomm) {
2341e07b27eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: preconditioner not yet setup\n");CHKERRQ(ierr);
2351e07b27eSBarry Smith     } else {
2361e07b27eSBarry Smith       MPI_Comm    comm,subcomm;
2371e07b27eSBarry Smith       PetscMPIInt comm_size,subcomm_size;
2381e07b27eSBarry Smith       DM          dm,subdm;
2391e07b27eSBarry Smith 
2401e07b27eSBarry Smith       ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
2411e07b27eSBarry Smith       subdm = private_PCTelescopeGetSubDM(sred);
2421e07b27eSBarry Smith       comm = PetscSubcommParent(sred->psubcomm);
2431e07b27eSBarry Smith       subcomm = PetscSubcommChild(sred->psubcomm);
2441e07b27eSBarry Smith       ierr = MPI_Comm_size(comm,&comm_size);CHKERRQ(ierr);
2451e07b27eSBarry Smith       ierr = MPI_Comm_size(subcomm,&subcomm_size);CHKERRQ(ierr);
2461e07b27eSBarry Smith 
2471e07b27eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: parent comm size reduction factor = %D\n",sred->redfactor);CHKERRQ(ierr);
2481e07b27eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: comm_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size);CHKERRQ(ierr);
249*48a10b22SPatrick Sanan       switch (sred->subcommtype) {
250*48a10b22SPatrick Sanan         case PETSC_SUBCOMM_INTERLACED :
251*48a10b22SPatrick Sanan           ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: subcomm type: interlaced\n",sred->subcommtype);CHKERRQ(ierr);
252*48a10b22SPatrick Sanan           break;
253*48a10b22SPatrick Sanan         case PETSC_SUBCOMM_CONTIGUOUS :
254*48a10b22SPatrick Sanan           ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: subcomm type: contiguous\n",sred->subcommtype);CHKERRQ(ierr);
255*48a10b22SPatrick Sanan           break;
256*48a10b22SPatrick Sanan         default :
257*48a10b22SPatrick Sanan           SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"General subcomm type not supported by PCTelescope");
258*48a10b22SPatrick Sanan       }
2591e07b27eSBarry Smith       ierr = PetscViewerGetSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
2601e07b27eSBarry Smith       if (isActiveRank(sred->psubcomm)) {
2611e07b27eSBarry Smith         ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr);
2621e07b27eSBarry Smith 
2631e07b27eSBarry Smith         if (dm && sred->ignore_dm) {
2641e07b27eSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: ignoring DM\n");CHKERRQ(ierr);
2651e07b27eSBarry Smith         }
2667c5279cbSDave May         if (sred->ignore_kspcomputeoperators) {
2677c5279cbSDave May           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: ignoring KSPComputeOperators\n");CHKERRQ(ierr);
2687c5279cbSDave May         }
2691e07b27eSBarry Smith         switch (sred->sr_type) {
2701e07b27eSBarry Smith         case TELESCOPE_DEFAULT:
2711e07b27eSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: using default setup\n");CHKERRQ(ierr);
2721e07b27eSBarry Smith           break;
2731e07b27eSBarry Smith         case TELESCOPE_DMDA:
2741e07b27eSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: DMDA detected\n");CHKERRQ(ierr);
2751e07b27eSBarry Smith           ierr = DMView_DMDAShort(subdm,subviewer);CHKERRQ(ierr);
2761e07b27eSBarry Smith           break;
2771e07b27eSBarry Smith         case TELESCOPE_DMPLEX:
2781e07b27eSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: DMPLEX detected\n");CHKERRQ(ierr);
2791e07b27eSBarry Smith           break;
2801e07b27eSBarry Smith         }
2811e07b27eSBarry Smith 
2821e07b27eSBarry Smith         ierr = KSPView(sred->ksp,subviewer);CHKERRQ(ierr);
2831e07b27eSBarry Smith         ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr);
2841e07b27eSBarry Smith       }
2851e07b27eSBarry Smith       ierr = PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
2861e07b27eSBarry Smith     }
2871e07b27eSBarry Smith   }
2881e07b27eSBarry Smith   PetscFunctionReturn(0);
2891e07b27eSBarry Smith }
2901e07b27eSBarry Smith 
2911e07b27eSBarry Smith #undef __FUNCT__
2921e07b27eSBarry Smith #define __FUNCT__ "PCSetUp_Telescope"
2931e07b27eSBarry Smith static PetscErrorCode PCSetUp_Telescope(PC pc)
2941e07b27eSBarry Smith {
2951e07b27eSBarry Smith   PC_Telescope      sred = (PC_Telescope)pc->data;
2961e07b27eSBarry Smith   PetscErrorCode    ierr;
297bd49479cSSatish Balay   MPI_Comm          comm,subcomm=0;
2981e07b27eSBarry Smith   PCTelescopeType   sr_type;
2991e07b27eSBarry Smith 
3001e07b27eSBarry Smith   PetscFunctionBegin;
3011e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
3021e07b27eSBarry Smith 
3031e07b27eSBarry Smith   /* subcomm definition */
3041e07b27eSBarry Smith   if (!pc->setupcalled) {
3051e07b27eSBarry Smith     if (!sred->psubcomm) {
3061e07b27eSBarry Smith       ierr = PetscSubcommCreate(comm,&sred->psubcomm);CHKERRQ(ierr);
3071e07b27eSBarry Smith       ierr = PetscSubcommSetNumber(sred->psubcomm,sred->redfactor);CHKERRQ(ierr);
308*48a10b22SPatrick Sanan       ierr = PetscSubcommSetType(sred->psubcomm,sred->subcommtype);CHKERRQ(ierr);
3091e07b27eSBarry Smith       ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr);
3101e07b27eSBarry Smith     }
3111e07b27eSBarry Smith   }
3120f6f40a7SSatish Balay   subcomm = PetscSubcommChild(sred->psubcomm);
3131e07b27eSBarry Smith 
3141e07b27eSBarry Smith   /* internal KSP */
3151e07b27eSBarry Smith   if (!pc->setupcalled) {
3161e07b27eSBarry Smith     const char *prefix;
3171e07b27eSBarry Smith 
3181e07b27eSBarry Smith     if (isActiveRank(sred->psubcomm)) {
3191e07b27eSBarry Smith       ierr = KSPCreate(subcomm,&sred->ksp);CHKERRQ(ierr);
3201e07b27eSBarry Smith       ierr = KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure);CHKERRQ(ierr);
3211e07b27eSBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
3221e07b27eSBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp);CHKERRQ(ierr);
3231e07b27eSBarry Smith       ierr = KSPSetType(sred->ksp,KSPPREONLY);CHKERRQ(ierr);
3241e07b27eSBarry Smith       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
3251e07b27eSBarry Smith       ierr = KSPSetOptionsPrefix(sred->ksp,prefix);CHKERRQ(ierr);
3261e07b27eSBarry Smith       ierr = KSPAppendOptionsPrefix(sred->ksp,"telescope_");CHKERRQ(ierr);
3271e07b27eSBarry Smith     }
3281e07b27eSBarry Smith   }
3291e07b27eSBarry Smith   /* Determine type of setup/update */
3301e07b27eSBarry Smith   if (!pc->setupcalled) {
3311e07b27eSBarry Smith     PetscBool has_dm,same;
3321e07b27eSBarry Smith     DM        dm;
3331e07b27eSBarry Smith 
3341e07b27eSBarry Smith     sr_type = TELESCOPE_DEFAULT;
3351e07b27eSBarry Smith     has_dm = PETSC_FALSE;
3361e07b27eSBarry Smith     ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
3371e07b27eSBarry Smith     if (dm) { has_dm = PETSC_TRUE; }
3381e07b27eSBarry Smith     if (has_dm) {
3391e07b27eSBarry Smith       /* check for dmda */
3401e07b27eSBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&same);CHKERRQ(ierr);
3411e07b27eSBarry Smith       if (same) {
3421e07b27eSBarry Smith         ierr = PetscInfo(pc,"PCTelescope: found DMDA\n");CHKERRQ(ierr);
3431e07b27eSBarry Smith         sr_type = TELESCOPE_DMDA;
3441e07b27eSBarry Smith       }
3451e07b27eSBarry Smith       /* check for dmplex */
3461e07b27eSBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same);CHKERRQ(ierr);
3471e07b27eSBarry Smith       if (same) {
3481e07b27eSBarry Smith         PetscInfo(pc,"PCTelescope: found DMPLEX\n");
3491e07b27eSBarry Smith         sr_type = TELESCOPE_DMPLEX;
3501e07b27eSBarry Smith       }
3511e07b27eSBarry Smith     }
3521e07b27eSBarry Smith 
3531e07b27eSBarry Smith     if (sred->ignore_dm) {
3541e07b27eSBarry Smith       PetscInfo(pc,"PCTelescope: ignore DM\n");
3551e07b27eSBarry Smith       sr_type = TELESCOPE_DEFAULT;
3561e07b27eSBarry Smith     }
3571e07b27eSBarry Smith     sred->sr_type = sr_type;
3581e07b27eSBarry Smith   } else {
3591e07b27eSBarry Smith     sr_type = sred->sr_type;
3601e07b27eSBarry Smith   }
3611e07b27eSBarry Smith 
3621e07b27eSBarry Smith   /* set function pointers for repartition setup, matrix creation/update, matrix nullspace and reset functionality */
3631e07b27eSBarry Smith   switch (sr_type) {
3641e07b27eSBarry Smith   case TELESCOPE_DEFAULT:
3651e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
3661e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
3671e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
3681e07b27eSBarry Smith     sred->pctelescope_reset_type              = NULL;
3691e07b27eSBarry Smith     break;
3701e07b27eSBarry Smith   case TELESCOPE_DMDA:
3711e07b27eSBarry Smith     pc->ops->apply                            = PCApply_Telescope_dmda;
372f650675bSDave May     pc->ops->applyrichardson                  = PCApplyRichardson_Telescope_dmda;
3731e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_dmda;
3741e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_dmda;
3751e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda;
3761e07b27eSBarry Smith     sred->pctelescope_reset_type              = PCReset_Telescope_dmda;
3771e07b27eSBarry Smith     break;
3781e07b27eSBarry Smith   case TELESCOPE_DMPLEX: SETERRQ(comm,PETSC_ERR_SUP,"Supprt for DMPLEX is currently not available");
3791e07b27eSBarry Smith     break;
3801e07b27eSBarry Smith   default: SETERRQ(comm,PETSC_ERR_SUP,"Only supprt for repartitioning DMDA is provided");
3811e07b27eSBarry Smith     break;
3821e07b27eSBarry Smith   }
3831e07b27eSBarry Smith 
3841e07b27eSBarry Smith   /* setup */
3851e07b27eSBarry Smith   if (sred->pctelescope_setup_type) {
3861e07b27eSBarry Smith     ierr = sred->pctelescope_setup_type(pc,sred);CHKERRQ(ierr);
3871e07b27eSBarry Smith   }
3881e07b27eSBarry Smith   /* update */
3891e07b27eSBarry Smith   if (!pc->setupcalled) {
3901e07b27eSBarry Smith     if (sred->pctelescope_matcreate_type) {
3911e07b27eSBarry Smith       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred);CHKERRQ(ierr);
3921e07b27eSBarry Smith     }
3931e07b27eSBarry Smith     if (sred->pctelescope_matnullspacecreate_type) {
3941e07b27eSBarry Smith       ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr);
3951e07b27eSBarry Smith     }
3961e07b27eSBarry Smith   } else {
3971e07b27eSBarry Smith     if (sred->pctelescope_matcreate_type) {
3981e07b27eSBarry Smith       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred);CHKERRQ(ierr);
3991e07b27eSBarry Smith     }
4001e07b27eSBarry Smith   }
4011e07b27eSBarry Smith 
4021e07b27eSBarry Smith   /* common - no construction */
4031e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
4041e07b27eSBarry Smith     ierr = KSPSetOperators(sred->ksp,sred->Bred,sred->Bred);CHKERRQ(ierr);
4051e07b27eSBarry Smith     if (pc->setfromoptionscalled && !pc->setupcalled){
4061e07b27eSBarry Smith       ierr = KSPSetFromOptions(sred->ksp);CHKERRQ(ierr);
4071e07b27eSBarry Smith     }
4081e07b27eSBarry Smith   }
4091e07b27eSBarry Smith   PetscFunctionReturn(0);
4101e07b27eSBarry Smith }
4111e07b27eSBarry Smith 
4121e07b27eSBarry Smith #undef __FUNCT__
4131e07b27eSBarry Smith #define __FUNCT__ "PCApply_Telescope"
4141e07b27eSBarry Smith static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y)
4151e07b27eSBarry Smith {
4161e07b27eSBarry Smith   PC_Telescope      sred = (PC_Telescope)pc->data;
4171e07b27eSBarry Smith   PetscErrorCode    ierr;
4181e07b27eSBarry Smith   Vec               xtmp,xred,yred;
41913c30530SDave May   PetscInt          i,st,ed;
4201e07b27eSBarry Smith   VecScatter        scatter;
4211e07b27eSBarry Smith   PetscScalar       *array;
4221e07b27eSBarry Smith   const PetscScalar *x_array;
4231e07b27eSBarry Smith 
4241e07b27eSBarry Smith   PetscFunctionBegin;
4251e07b27eSBarry Smith   xtmp    = sred->xtmp;
4261e07b27eSBarry Smith   scatter = sred->scatter;
4271e07b27eSBarry Smith   xred    = sred->xred;
4281e07b27eSBarry Smith   yred    = sred->yred;
4291e07b27eSBarry Smith 
4301e07b27eSBarry Smith   /* pull in vector x->xtmp */
4311e07b27eSBarry Smith   ierr = VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4321e07b27eSBarry Smith   ierr = VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4331e07b27eSBarry Smith 
4341e07b27eSBarry Smith   /* copy vector entires into xred */
4351e07b27eSBarry Smith   ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr);
4361e07b27eSBarry Smith   if (xred) {
4371e07b27eSBarry Smith     PetscScalar *LA_xred;
4381e07b27eSBarry Smith     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
4391e07b27eSBarry Smith     ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr);
4401e07b27eSBarry Smith     for (i=0; i<ed-st; i++) {
4411e07b27eSBarry Smith       LA_xred[i] = x_array[i];
4421e07b27eSBarry Smith     }
4431e07b27eSBarry Smith     ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr);
4441e07b27eSBarry Smith   }
4451e07b27eSBarry Smith   ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr);
4461e07b27eSBarry Smith   /* solve */
4471e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
4481e07b27eSBarry Smith     ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr);
4491e07b27eSBarry Smith   }
4501e07b27eSBarry Smith   /* return vector */
4511e07b27eSBarry Smith   ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr);
4521e07b27eSBarry Smith   if (yred) {
4531e07b27eSBarry Smith     const PetscScalar *LA_yred;
4541e07b27eSBarry Smith     ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr);
4551e07b27eSBarry Smith     ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr);
4561e07b27eSBarry Smith     for (i=0; i<ed-st; i++) {
4571e07b27eSBarry Smith       array[i] = LA_yred[i];
4581e07b27eSBarry Smith     }
4591e07b27eSBarry Smith     ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr);
4601e07b27eSBarry Smith   }
4611e07b27eSBarry Smith   ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr);
4621e07b27eSBarry Smith   ierr = VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4631e07b27eSBarry Smith   ierr = VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4641e07b27eSBarry Smith   PetscFunctionReturn(0);
4651e07b27eSBarry Smith }
4661e07b27eSBarry Smith 
4671e07b27eSBarry Smith #undef __FUNCT__
468f650675bSDave May #define __FUNCT__ "PCApplyRichardson_Telescope"
469f650675bSDave 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)
470f650675bSDave May {
471f650675bSDave May   PC_Telescope      sred = (PC_Telescope)pc->data;
472f650675bSDave May   PetscErrorCode    ierr;
473a1d91a28SDave May   Vec               xtmp,yred;
474f650675bSDave May   PetscInt          i,st,ed;
475f650675bSDave May   VecScatter        scatter;
476f650675bSDave May   const PetscScalar *x_array;
477f650675bSDave May   PetscBool         default_init_guess_value;
478f650675bSDave May 
479f650675bSDave May   PetscFunctionBegin;
480f650675bSDave May   xtmp    = sred->xtmp;
481f650675bSDave May   scatter = sred->scatter;
482f650675bSDave May   yred    = sred->yred;
483f650675bSDave May 
484f650675bSDave May   if (its > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCApplyRichardson_Telescope only supports max_it = 1");
485f650675bSDave May   *reason = (PCRichardsonConvergedReason)0;
486f650675bSDave May 
487f650675bSDave May   if (!zeroguess) {
488f650675bSDave May     ierr = PetscInfo(pc,"PCTelescope: Scattering y for non-zero initial guess\n");CHKERRQ(ierr);
489f650675bSDave May     /* pull in vector y->xtmp */
490f650675bSDave May     ierr = VecScatterBegin(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
491f650675bSDave May     ierr = VecScatterEnd(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
492f650675bSDave May 
493f650675bSDave May     /* copy vector entires into xred */
494f650675bSDave May     ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr);
495f650675bSDave May     if (yred) {
496f650675bSDave May       PetscScalar *LA_yred;
497f650675bSDave May       ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr);
498f650675bSDave May       ierr = VecGetArray(yred,&LA_yred);CHKERRQ(ierr);
499f650675bSDave May       for (i=0; i<ed-st; i++) {
500f650675bSDave May         LA_yred[i] = x_array[i];
501f650675bSDave May       }
502f650675bSDave May       ierr = VecRestoreArray(yred,&LA_yred);CHKERRQ(ierr);
503f650675bSDave May     }
504f650675bSDave May     ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr);
505f650675bSDave May   }
506f650675bSDave May 
507f650675bSDave May   if (isActiveRank(sred->psubcomm)) {
508f650675bSDave May     ierr = KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value);CHKERRQ(ierr);
509f650675bSDave May     if (!zeroguess) ierr = KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE);CHKERRQ(ierr);
510f650675bSDave May   }
511f650675bSDave May 
512f650675bSDave May   ierr = PCApply_Telescope(pc,x,y);CHKERRQ(ierr);
513f650675bSDave May 
514f650675bSDave May   if (isActiveRank(sred->psubcomm)) {
515f650675bSDave May     ierr = KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value);CHKERRQ(ierr);
516f650675bSDave May   }
517f650675bSDave May 
518f650675bSDave May   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
519f650675bSDave May   *outits = 1;
520f650675bSDave May   PetscFunctionReturn(0);
521f650675bSDave May }
522f650675bSDave May 
523f650675bSDave May #undef __FUNCT__
5241e07b27eSBarry Smith #define __FUNCT__ "PCReset_Telescope"
5251e07b27eSBarry Smith static PetscErrorCode PCReset_Telescope(PC pc)
5261e07b27eSBarry Smith {
5271e07b27eSBarry Smith   PC_Telescope   sred = (PC_Telescope)pc->data;
5281e07b27eSBarry Smith   PetscErrorCode ierr;
5291e07b27eSBarry Smith 
5301e07b27eSBarry Smith   ierr = ISDestroy(&sred->isin);CHKERRQ(ierr);
5311e07b27eSBarry Smith   ierr = VecScatterDestroy(&sred->scatter);CHKERRQ(ierr);
532e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->xred);CHKERRQ(ierr);
533e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->yred);CHKERRQ(ierr);
534e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->xtmp);CHKERRQ(ierr);
535e3acf2f7SBarry Smith   ierr = MatDestroy(&sred->Bred);CHKERRQ(ierr);
536e3acf2f7SBarry Smith   ierr = KSPReset(sred->ksp);CHKERRQ(ierr);
5371e07b27eSBarry Smith   if (sred->pctelescope_reset_type) {
5381e07b27eSBarry Smith     ierr = sred->pctelescope_reset_type(pc);CHKERRQ(ierr);
5391e07b27eSBarry Smith   }
5401e07b27eSBarry Smith   PetscFunctionReturn(0);
5411e07b27eSBarry Smith }
5421e07b27eSBarry Smith 
5431e07b27eSBarry Smith #undef __FUNCT__
5441e07b27eSBarry Smith #define __FUNCT__ "PCDestroy_Telescope"
5451e07b27eSBarry Smith static PetscErrorCode PCDestroy_Telescope(PC pc)
5461e07b27eSBarry Smith {
5471e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
5481e07b27eSBarry Smith   PetscErrorCode   ierr;
5491e07b27eSBarry Smith 
5501e07b27eSBarry Smith   PetscFunctionBegin;
5511e07b27eSBarry Smith   ierr = PCReset_Telescope(pc);CHKERRQ(ierr);
552e3acf2f7SBarry Smith   ierr = KSPDestroy(&sred->ksp);CHKERRQ(ierr);
5531e07b27eSBarry Smith   ierr = PetscSubcommDestroy(&sred->psubcomm);CHKERRQ(ierr);
554e3acf2f7SBarry Smith   ierr = PetscFree(sred->dm_ctx);CHKERRQ(ierr);
555e3acf2f7SBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
5561e07b27eSBarry Smith   PetscFunctionReturn(0);
5571e07b27eSBarry Smith }
5581e07b27eSBarry Smith 
5591e07b27eSBarry Smith #undef __FUNCT__
5601e07b27eSBarry Smith #define __FUNCT__ "PCSetFromOptions_Telescope"
5614416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Telescope(PetscOptionItems *PetscOptionsObject,PC pc)
5621e07b27eSBarry Smith {
5631e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
5641e07b27eSBarry Smith   PetscErrorCode   ierr;
5651e07b27eSBarry Smith   MPI_Comm         comm;
5661e07b27eSBarry Smith   PetscMPIInt      size;
567*48a10b22SPatrick Sanan   PetscBool        flg;
568*48a10b22SPatrick Sanan   PetscSubcommType subcommtype;
5691e07b27eSBarry Smith 
5701e07b27eSBarry Smith   PetscFunctionBegin;
5711e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
5721e07b27eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
5731e07b27eSBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Telescope options");CHKERRQ(ierr);
574*48a10b22SPatrick Sanan   ierr = PetscOptionsEnum("-pc_telescope_subcomm_type","Subcomm type (interlaced or contiguous)","PCTelescopeSetSubcommType",PetscSubcommTypes,(PetscEnum)sred->subcommtype,(PetscEnum*)&subcommtype,&flg);CHKERRQ(ierr);
575*48a10b22SPatrick Sanan   if (flg) {
576*48a10b22SPatrick Sanan     ierr = PCTelescopeSetSubcommType(pc,subcommtype);CHKERRQ(ierr);
577*48a10b22SPatrick Sanan   }
5781e07b27eSBarry Smith   ierr = PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,0);CHKERRQ(ierr);
5791e07b27eSBarry Smith   if (sred->redfactor > size) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size");
5801e07b27eSBarry Smith   ierr = PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,0);CHKERRQ(ierr);
5817c5279cbSDave May   ierr = PetscOptionsBool("-pc_telescope_ignore_kspcomputeoperators","Ignore method used to compute A","PCTelescopeSetIgnoreKSPComputeOperators",sred->ignore_kspcomputeoperators,&sred->ignore_kspcomputeoperators,0);CHKERRQ(ierr);
5821e07b27eSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
5831e07b27eSBarry Smith   PetscFunctionReturn(0);
5841e07b27eSBarry Smith }
5851e07b27eSBarry Smith 
5861e07b27eSBarry Smith /* PC simplementation specific API's */
5871e07b27eSBarry Smith 
588*48a10b22SPatrick Sanan #undef __FUNCT__
589*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetKSP_Telescope"
5901e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp)
5911e07b27eSBarry Smith {
5921e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
593bd49479cSSatish Balay   PetscFunctionBegin;
5941e07b27eSBarry Smith   if (ksp) *ksp = red->ksp;
595bd49479cSSatish Balay   PetscFunctionReturn(0);
5961e07b27eSBarry Smith }
5971e07b27eSBarry Smith 
598*48a10b22SPatrick Sanan #undef __FUNCT__
599*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetSubcommType_Telescope"
600*48a10b22SPatrick Sanan static PetscErrorCode PCTelescopeGetSubcommType_Telescope(PC pc,PetscSubcommType *subcommtype)
601*48a10b22SPatrick Sanan {
602*48a10b22SPatrick Sanan   PC_Telescope red = (PC_Telescope)pc->data;
603*48a10b22SPatrick Sanan   PetscFunctionBegin;
604*48a10b22SPatrick Sanan   if (subcommtype) *subcommtype = red->subcommtype;
605*48a10b22SPatrick Sanan   PetscFunctionReturn(0);
606*48a10b22SPatrick Sanan }
607*48a10b22SPatrick Sanan 
608*48a10b22SPatrick Sanan #undef __FUNCT__
609*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeSetSubcommType_Telescope"
610*48a10b22SPatrick Sanan static PetscErrorCode PCTelescopeSetSubcommType_Telescope(PC pc,PetscSubcommType subcommtype)
611*48a10b22SPatrick Sanan {
612*48a10b22SPatrick Sanan   PC_Telescope     red = (PC_Telescope)pc->data;
613*48a10b22SPatrick Sanan 
614*48a10b22SPatrick Sanan   PetscFunctionBegin;
615*48a10b22SPatrick 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.");
616*48a10b22SPatrick Sanan   red->subcommtype = subcommtype;
617*48a10b22SPatrick Sanan   PetscFunctionReturn(0);
618*48a10b22SPatrick Sanan }
619*48a10b22SPatrick Sanan 
620*48a10b22SPatrick Sanan #undef __FUNCT__
621*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetReductionFactor_Telescope"
6221e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact)
6231e07b27eSBarry Smith {
6241e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
625bd49479cSSatish Balay   PetscFunctionBegin;
6261e07b27eSBarry Smith   if (fact) *fact = red->redfactor;
627bd49479cSSatish Balay   PetscFunctionReturn(0);
6281e07b27eSBarry Smith }
6291e07b27eSBarry Smith 
630*48a10b22SPatrick Sanan #undef __FUNCT__
631*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeSetReductionFactor_Telescope"
6321e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact)
6331e07b27eSBarry Smith {
6341e07b27eSBarry Smith   PC_Telescope     red = (PC_Telescope)pc->data;
6351e07b27eSBarry Smith   PetscMPIInt      size;
6361e07b27eSBarry Smith   PetscErrorCode   ierr;
6371e07b27eSBarry Smith 
638bd49479cSSatish Balay   PetscFunctionBegin;
6391e07b27eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
6401e07b27eSBarry Smith   if (fact <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact);
6411e07b27eSBarry Smith   if (fact > size) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact);
6421e07b27eSBarry Smith   red->redfactor = fact;
643bd49479cSSatish Balay   PetscFunctionReturn(0);
6441e07b27eSBarry Smith }
6451e07b27eSBarry Smith 
646*48a10b22SPatrick Sanan #undef __FUNCT__
647*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetIgnoreDM_Telescope"
6481e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v)
6491e07b27eSBarry Smith {
6501e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
651bd49479cSSatish Balay   PetscFunctionBegin;
6521e07b27eSBarry Smith   if (v) *v = red->ignore_dm;
653bd49479cSSatish Balay   PetscFunctionReturn(0);
6541e07b27eSBarry Smith }
655*48a10b22SPatrick Sanan 
656*48a10b22SPatrick Sanan #undef __FUNCT__
657*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeSetIgnoreDM_Telescope"
6581e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v)
6591e07b27eSBarry Smith {
6601e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
661bd49479cSSatish Balay   PetscFunctionBegin;
6621e07b27eSBarry Smith   red->ignore_dm = v;
663bd49479cSSatish Balay   PetscFunctionReturn(0);
6641e07b27eSBarry Smith }
6651e07b27eSBarry Smith 
666*48a10b22SPatrick Sanan #undef __FUNCT__
667*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetIgnoreKSPComputeOperators_Telescope"
6680ae7c45bSDave May static PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool *v)
6690ae7c45bSDave May {
6700ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
6710ae7c45bSDave May   PetscFunctionBegin;
6720ae7c45bSDave May   if (v) *v = red->ignore_kspcomputeoperators;
6730ae7c45bSDave May   PetscFunctionReturn(0);
6740ae7c45bSDave May }
675*48a10b22SPatrick Sanan 
676*48a10b22SPatrick Sanan #undef __FUNCT__
677*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeSetIgnoreKSPComputeOperators_Telescope"
6780ae7c45bSDave May static PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool v)
6790ae7c45bSDave May {
6800ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
6810ae7c45bSDave May   PetscFunctionBegin;
6820ae7c45bSDave May   red->ignore_kspcomputeoperators = v;
6830ae7c45bSDave May   PetscFunctionReturn(0);
6840ae7c45bSDave May }
6850ae7c45bSDave May 
686*48a10b22SPatrick Sanan #undef __FUNCT__
687*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetDM_Telescope"
6881e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm)
6891e07b27eSBarry Smith {
6901e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
691bd49479cSSatish Balay   PetscFunctionBegin;
6921e07b27eSBarry Smith   *dm = private_PCTelescopeGetSubDM(red);
693bd49479cSSatish Balay   PetscFunctionReturn(0);
6941e07b27eSBarry Smith }
6951e07b27eSBarry Smith 
6961e07b27eSBarry Smith /*@
6971e07b27eSBarry Smith  PCTelescopeGetKSP - Gets the KSP created by the telescoping PC.
6981e07b27eSBarry Smith 
6991e07b27eSBarry Smith  Not Collective
7001e07b27eSBarry Smith 
7011e07b27eSBarry Smith  Input Parameter:
7021e07b27eSBarry Smith  .  pc - the preconditioner context
7031e07b27eSBarry Smith 
7041e07b27eSBarry Smith  Output Parameter:
7051e07b27eSBarry Smith  .  subksp - the KSP defined the smaller set of processes
7061e07b27eSBarry Smith 
7071e07b27eSBarry Smith  Level: advanced
7081e07b27eSBarry Smith 
7091e07b27eSBarry Smith  .keywords: PC, telescopting solve
7101e07b27eSBarry Smith  @*/
711*48a10b22SPatrick Sanan #undef __FUNCT__
712*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetKSP"
7131e07b27eSBarry Smith PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp)
7141e07b27eSBarry Smith {
715bd49479cSSatish Balay   PetscErrorCode ierr;
716bd49479cSSatish Balay   PetscFunctionBegin;
717163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp));CHKERRQ(ierr);
718bd49479cSSatish Balay   PetscFunctionReturn(0);
7191e07b27eSBarry Smith }
7201e07b27eSBarry Smith 
7211e07b27eSBarry Smith /*@
7221e07b27eSBarry Smith  PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by.
7231e07b27eSBarry Smith 
7241e07b27eSBarry Smith  Not Collective
7251e07b27eSBarry Smith 
7261e07b27eSBarry Smith  Input Parameter:
7271e07b27eSBarry Smith  .  pc - the preconditioner context
7281e07b27eSBarry Smith 
7291e07b27eSBarry Smith  Output Parameter:
7301e07b27eSBarry Smith  .  fact - the reduction factor
7311e07b27eSBarry Smith 
7321e07b27eSBarry Smith  Level: advanced
7331e07b27eSBarry Smith 
7341e07b27eSBarry Smith  .keywords: PC, telescoping solve
7351e07b27eSBarry Smith  @*/
736*48a10b22SPatrick Sanan #undef __FUNCT__
737*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetReductionFactor"
7381e07b27eSBarry Smith PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact)
7391e07b27eSBarry Smith {
740bd49479cSSatish Balay   PetscErrorCode ierr;
741bd49479cSSatish Balay   PetscFunctionBegin;
742163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact));CHKERRQ(ierr);
743bd49479cSSatish Balay   PetscFunctionReturn(0);
7441e07b27eSBarry Smith }
7451e07b27eSBarry Smith 
7461e07b27eSBarry Smith /*@
7471e07b27eSBarry Smith  PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by.
7481e07b27eSBarry Smith 
7491e07b27eSBarry Smith  Not Collective
7501e07b27eSBarry Smith 
7511e07b27eSBarry Smith  Input Parameter:
7521e07b27eSBarry Smith  .  pc - the preconditioner context
7531e07b27eSBarry Smith 
7541e07b27eSBarry Smith  Output Parameter:
7551e07b27eSBarry Smith  .  fact - the reduction factor
7561e07b27eSBarry Smith 
7571e07b27eSBarry Smith  Level: advanced
7581e07b27eSBarry Smith 
7591e07b27eSBarry Smith  .keywords: PC, telescoping solve
7601e07b27eSBarry Smith  @*/
761*48a10b22SPatrick Sanan #undef __FUNCT__
762*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeSetReductionFactor"
7631e07b27eSBarry Smith PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact)
7641e07b27eSBarry Smith {
765bd49479cSSatish Balay   PetscErrorCode ierr;
766bd49479cSSatish Balay   PetscFunctionBegin;
767bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact));CHKERRQ(ierr);
768bd49479cSSatish Balay   PetscFunctionReturn(0);
7691e07b27eSBarry Smith }
7701e07b27eSBarry Smith 
7711e07b27eSBarry Smith /*@
7721e07b27eSBarry Smith  PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used.
7731e07b27eSBarry Smith 
7741e07b27eSBarry Smith  Not Collective
7751e07b27eSBarry Smith 
7761e07b27eSBarry Smith  Input Parameter:
7771e07b27eSBarry Smith  .  pc - the preconditioner context
7781e07b27eSBarry Smith 
7791e07b27eSBarry Smith  Output Parameter:
7801e07b27eSBarry Smith  .  v - the flag
7811e07b27eSBarry Smith 
7821e07b27eSBarry Smith  Level: advanced
7831e07b27eSBarry Smith 
7841e07b27eSBarry Smith  .keywords: PC, telescoping solve
7851e07b27eSBarry Smith  @*/
786*48a10b22SPatrick Sanan #undef __FUNCT__
787*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetIgnoreDM"
7881e07b27eSBarry Smith PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v)
7891e07b27eSBarry Smith {
790bd49479cSSatish Balay   PetscErrorCode ierr;
791bd49479cSSatish Balay   PetscFunctionBegin;
792163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr);
793bd49479cSSatish Balay   PetscFunctionReturn(0);
7941e07b27eSBarry Smith }
7951e07b27eSBarry Smith 
7961e07b27eSBarry Smith /*@
7971e07b27eSBarry Smith  PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC.
7981e07b27eSBarry Smith 
7991e07b27eSBarry Smith  Not Collective
8001e07b27eSBarry Smith 
8011e07b27eSBarry Smith  Input Parameter:
8021e07b27eSBarry Smith  .  pc - the preconditioner context
8031e07b27eSBarry Smith 
8041e07b27eSBarry Smith  Output Parameter:
8051e07b27eSBarry Smith  .  v - Use PETSC_TRUE to ignore any DM
8061e07b27eSBarry Smith 
8071e07b27eSBarry Smith  Level: advanced
8081e07b27eSBarry Smith 
8091e07b27eSBarry Smith  .keywords: PC, telescoping solve
8101e07b27eSBarry Smith  @*/
811*48a10b22SPatrick Sanan #undef __FUNCT__
812*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeSetIgnoreDM"
813bfd6bcc6SSatish Balay PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v)
8141e07b27eSBarry Smith {
815bd49479cSSatish Balay   PetscErrorCode ierr;
816bd49479cSSatish Balay   PetscFunctionBegin;
817bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr);
818bd49479cSSatish Balay   PetscFunctionReturn(0);
8191e07b27eSBarry Smith }
8201e07b27eSBarry Smith 
8211e07b27eSBarry Smith /*@
8220ae7c45bSDave May  PCTelescopeGetIgnoreKSPComputeOperators - Get the flag indicating if KSPComputeOperators will be used.
8230ae7c45bSDave May 
8240ae7c45bSDave May  Not Collective
8250ae7c45bSDave May 
8260ae7c45bSDave May  Input Parameter:
8270ae7c45bSDave May  .  pc - the preconditioner context
8280ae7c45bSDave May 
8290ae7c45bSDave May  Output Parameter:
8300ae7c45bSDave May  .  v - the flag
8310ae7c45bSDave May 
8320ae7c45bSDave May  Level: advanced
8330ae7c45bSDave May 
8340ae7c45bSDave May  .keywords: PC, telescoping solve
8350ae7c45bSDave May  @*/
836*48a10b22SPatrick Sanan #undef __FUNCT__
837*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetIgnoreKSPComputeOperators"
8380ae7c45bSDave May PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC pc,PetscBool *v)
8390ae7c45bSDave May {
8400ae7c45bSDave May   PetscErrorCode ierr;
8410ae7c45bSDave May   PetscFunctionBegin;
842163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr);
8430ae7c45bSDave May   PetscFunctionReturn(0);
8440ae7c45bSDave May }
8450ae7c45bSDave May 
8460ae7c45bSDave May /*@
8470ae7c45bSDave May  PCTelescopeSetIgnoreKSPComputeOperators - Set a flag to ignore KSPComputeOperators.
8480ae7c45bSDave May 
8490ae7c45bSDave May  Not Collective
8500ae7c45bSDave May 
8510ae7c45bSDave May  Input Parameter:
8520ae7c45bSDave May  .  pc - the preconditioner context
8530ae7c45bSDave May 
8540ae7c45bSDave May  Output Parameter:
855a954d8f4SDave May  .  v - Use PETSC_TRUE to ignore the method (if defined) set via KSPSetComputeOperators on pc
8560ae7c45bSDave May 
8570ae7c45bSDave May  Level: advanced
8580ae7c45bSDave May 
8590ae7c45bSDave May  .keywords: PC, telescoping solve
8600ae7c45bSDave May  @*/
861*48a10b22SPatrick Sanan #undef __FUNCT__
862*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeSetIgnoreKSPComputeOperators"
8630ae7c45bSDave May PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC pc,PetscBool v)
8640ae7c45bSDave May {
8650ae7c45bSDave May   PetscErrorCode ierr;
8660ae7c45bSDave May   PetscFunctionBegin;
8670ae7c45bSDave May   ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr);
8680ae7c45bSDave May   PetscFunctionReturn(0);
8690ae7c45bSDave May }
8700ae7c45bSDave May 
8710ae7c45bSDave May /*@
8721e07b27eSBarry Smith  PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP.
8731e07b27eSBarry Smith 
8741e07b27eSBarry Smith  Not Collective
8751e07b27eSBarry Smith 
8761e07b27eSBarry Smith  Input Parameter:
8771e07b27eSBarry Smith  .  pc - the preconditioner context
8781e07b27eSBarry Smith 
8791e07b27eSBarry Smith  Output Parameter:
8801e07b27eSBarry Smith  .  subdm - The re-partitioned DM
8811e07b27eSBarry Smith 
8821e07b27eSBarry Smith  Level: advanced
8831e07b27eSBarry Smith 
8841e07b27eSBarry Smith  .keywords: PC, telescoping solve
8851e07b27eSBarry Smith  @*/
886*48a10b22SPatrick Sanan #undef __FUNCT__
887*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetDM"
8881e07b27eSBarry Smith PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm)
8891e07b27eSBarry Smith {
890bd49479cSSatish Balay   PetscErrorCode ierr;
891bd49479cSSatish Balay   PetscFunctionBegin;
892163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm));CHKERRQ(ierr);
893bd49479cSSatish Balay   PetscFunctionReturn(0);
8941e07b27eSBarry Smith }
8951e07b27eSBarry Smith 
896*48a10b22SPatrick Sanan /*@
897*48a10b22SPatrick Sanan  PCTelescopeSetSubcommType - set subcommunicator type (interlaced or contiguous)
898*48a10b22SPatrick Sanan 
899*48a10b22SPatrick Sanan  Logically Collective
900*48a10b22SPatrick Sanan 
901*48a10b22SPatrick Sanan  Input Parameter:
902*48a10b22SPatrick Sanan  .  pc - the preconditioner context
903*48a10b22SPatrick Sanan  .  subcommtype - the subcommunicator type (see PetscSubcommType)
904*48a10b22SPatrick Sanan 
905*48a10b22SPatrick Sanan  Level: advanced
906*48a10b22SPatrick Sanan 
907*48a10b22SPatrick Sanan  .keywords: PC, telescoping solve
908*48a10b22SPatrick Sanan 
909*48a10b22SPatrick Sanan  .seealso : PetscSubcommType, PetscSubcomm, PCTELESCOPE
910*48a10b22SPatrick Sanan  @*/
911*48a10b22SPatrick Sanan #undef __FUNCT__
912*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeSetSubcommType"
913*48a10b22SPatrick Sanan PetscErrorCode PCTelescopeSetSubcommType(PC pc, PetscSubcommType subcommtype)
914*48a10b22SPatrick Sanan {
915*48a10b22SPatrick Sanan   PetscErrorCode ierr;
916*48a10b22SPatrick Sanan   PetscFunctionBegin;
917*48a10b22SPatrick Sanan   ierr = PetscTryMethod(pc,"PCTelescopeSetSubcommType_C",(PC,PetscSubcommType),(pc,subcommtype));CHKERRQ(ierr);
918*48a10b22SPatrick Sanan   PetscFunctionReturn(0);
919*48a10b22SPatrick Sanan }
920*48a10b22SPatrick Sanan 
921*48a10b22SPatrick Sanan /*@
922*48a10b22SPatrick Sanan  PCTelescopeGetSubcommType - Get the subcommunicator type (interlaced or contiguous)
923*48a10b22SPatrick Sanan 
924*48a10b22SPatrick Sanan  Not Collective
925*48a10b22SPatrick Sanan 
926*48a10b22SPatrick Sanan  Input Parameter:
927*48a10b22SPatrick Sanan  .  pc - the preconditioner context
928*48a10b22SPatrick Sanan 
929*48a10b22SPatrick Sanan  Output Parameter:
930*48a10b22SPatrick Sanan  .  subcommtype - the subcommunicator type (see PetscSubcommType)
931*48a10b22SPatrick Sanan 
932*48a10b22SPatrick Sanan  Level: advanced
933*48a10b22SPatrick Sanan 
934*48a10b22SPatrick Sanan  .keywords: PC, telescoping solve
935*48a10b22SPatrick Sanan 
936*48a10b22SPatrick Sanan  .seealso : PetscSubComm, PetscSubcommType, PCTELESCOPE
937*48a10b22SPatrick Sanan  @*/
938*48a10b22SPatrick Sanan #undef __FUNCT__
939*48a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetSubCommType"
940*48a10b22SPatrick Sanan PetscErrorCode PCTelescopeGetSubCommType(PC pc, PetscSubcommType *subcommtype)
941*48a10b22SPatrick Sanan {
942*48a10b22SPatrick Sanan   PetscErrorCode ierr;
943*48a10b22SPatrick Sanan   PetscFunctionBegin;
944*48a10b22SPatrick Sanan   ierr = PetscUseMethod(pc,"PCTelescopeGetSubcommType_C",(PC,PetscSubcommType*),(pc,subcommtype));CHKERRQ(ierr);
945*48a10b22SPatrick Sanan   PetscFunctionReturn(0);
946*48a10b22SPatrick Sanan }
947*48a10b22SPatrick Sanan 
9481e07b27eSBarry Smith /* -------------------------------------------------------------------------------------*/
9491e07b27eSBarry Smith /*MC
9501e07b27eSBarry Smith    PCTELESCOPE - Runs a KSP solver on a sub-group of processors. MPI processes not in the sub-communicator are idle during the solve.
9511e07b27eSBarry Smith 
9521e07b27eSBarry Smith    Options Database:
9531e07b27eSBarry Smith +  -pc_telescope_reduction_factor <n> - factor to use communicator size by, for example if you are using 64 MPI processes and
9541e07b27eSBarry Smith    use an n of 4, the new sub-communicator will be 4 defined with 64/4 processes
9551e07b27eSBarry Smith -  -pc_telescope_ignore_dm <false> - flag to indicate whether an attached DM should be ignored
9561e07b27eSBarry Smith 
9571e07b27eSBarry Smith    Level: advanced
9581e07b27eSBarry Smith 
9591e07b27eSBarry Smith    Notes:
9606fc41876SBarry Smith    The preconditioner is deemed telescopic as it only calls KSPSolve() on a single
9616fc41876SBarry Smith    sub-communicator in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators.
9626fc41876SBarry Smith    This means there will be MPI processes within c, which will be idle during the application of this preconditioner.
9636fc41876SBarry Smith 
9641e07b27eSBarry Smith    The default KSP is PREONLY. If a DM is attached to the PC, it is re-partitioned on the sub-communicator.
9651e07b27eSBarry Smith    Both the B mat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned DM.
9661e07b27eSBarry Smith    Currently only support for re-partitioning a DMDA is provided.
9671e07b27eSBarry Smith    Any nullspace attached to the original Bmat operator are extracted, re-partitioned and set on the repartitioned Bmat operator.
9681e07b27eSBarry Smith    KSPSetComputeOperators() is not propagated to the sub KSP.
9691e07b27eSBarry Smith    Currently there is no support for the flag -pc_use_amat
9701e07b27eSBarry Smith 
9716fc41876SBarry Smith    Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation
9726fc41876SBarry Smith    creates a child sub-communicator (c') containing less MPI processes than the original parent preconditioner (PC).
9736fc41876SBarry Smith 
9746fc41876SBarry Smith   Developer Notes:
9756fc41876SBarry Smith    During PCSetup, the B operator is scattered onto c'.
9766fc41876SBarry Smith    Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c').
9776fc41876SBarry Smith    Then KSPSolve() is executed on the c' communicator.
9786fc41876SBarry Smith 
9796fc41876SBarry Smith    The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED
9806fc41876SBarry Smith    creation routine. We run the sub KSP on only the ranks within the communicator which have a color equal to zero.
9816fc41876SBarry Smith 
9826fc41876SBarry Smith    The telescoping preconditioner is aware of nullspaces which are attached to the only B operator.
9836fc41876SBarry Smith    In case where B has a n nullspace attached, these nullspaces vectors are extract from B and mapped into
9846fc41876SBarry Smith    a new nullspace (defined on the sub-communicator) which is attached to B' (the B operator which was scattered to c')
9856fc41876SBarry Smith 
9866fc41876SBarry Smith    The telescoping preconditioner is aware of an attached DM. In the event that the DM is of type DMDA (2D or 3D -
9876fc41876SBarry 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
9886fc41876SBarry Smith    is attached the sub KSPSolve(). The design of telescope is such that it should be possible to extend support
9896fc41876SBarry Smith    for re-partitioning other DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs.
9906fc41876SBarry Smith 
9916fc41876SBarry Smith    By default, B' is defined by simply fusing rows from different MPI processes
9926fc41876SBarry Smith 
9936fc41876SBarry Smith    When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permuting of B
9946fc41876SBarry Smith    into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatGetSubMatrices(), (iii) fusing the
9956fc41876SBarry Smith    locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat()
9966fc41876SBarry Smith 
9976fc41876SBarry Smith    Limitations/improvements
9986fc41876SBarry Smith    VecPlaceArray could be used within PCApply() to improve efficiency and reduce memory usage.
9996fc41876SBarry Smith 
10006fc41876SBarry Smith    The symmetric permutation used when a DMDA is encountered is performed via explicitly assmbleming a permutation matrix P,
10016fc41876SBarry Smith    and performing P^T.A.P. Possibly it might be more efficient to use MatPermute(). I opted to use P^T.A.P as it appears
10026fc41876SBarry Smith    VecPermute() does not supported for the use case required here. By computing P, I can permute both the operator and RHS in a
10036fc41876SBarry Smith    consistent manner.
10046fc41876SBarry Smith 
10056fc41876SBarry Smith    Mapping of vectors is performed this way
10066fc41876SBarry Smith    Suppose the parent comm size was 4, and we set a reduction factor of 2, thus would give a comm size on c' of 2.
10076fc41876SBarry Smith    Using the interlaced creation routine, the ranks in c with color = 0, will be rank 0 and 2.
10086fc41876SBarry Smith    We perform the scatter to the sub-comm in the following way,
10096fc41876SBarry Smith    [1] Given a vector x defined on comm c
10106fc41876SBarry Smith 
10116fc41876SBarry Smith    rank(c) : _________ 0 ______  ________ 1 _______  ________ 2 _____________ ___________ 3 __________
10126fc41876SBarry 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]
10136fc41876SBarry Smith 
10146fc41876SBarry Smith    scatter to xtmp defined also on comm c so that we have the following values
10156fc41876SBarry Smith 
10166fc41876SBarry Smith    rank(c) : ___________________ 0 ________________  _1_  ______________________ 2 _______________________  __3_
10176fc41876SBarry 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] [  ]
10186fc41876SBarry Smith 
10196fc41876SBarry Smith    The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') have no values
10206fc41876SBarry Smith 
10216fc41876SBarry Smith 
10226fc41876SBarry 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'.
10236fc41876SBarry Smith    Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0.
10246fc41876SBarry Smith 
10256fc41876SBarry Smith     rank(c') : ___________________ 0 _______________  ______________________ 1 _____________________
10266fc41876SBarry 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]
10276fc41876SBarry Smith 
10281e07b27eSBarry Smith 
10291e07b27eSBarry Smith   Contributed by Dave May
10301e07b27eSBarry Smith 
10316fc41876SBarry Smith .seealso:  PCTelescopeGetKSP(), PCTelescopeGetDM(), PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(), PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT
10321e07b27eSBarry Smith M*/
10331e07b27eSBarry Smith #undef __FUNCT__
10341e07b27eSBarry Smith #define __FUNCT__ "PCCreate_Telescope"
10351e07b27eSBarry Smith PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc)
10361e07b27eSBarry Smith {
10371e07b27eSBarry Smith   PetscErrorCode       ierr;
10381e07b27eSBarry Smith   struct _PC_Telescope *sred;
10391e07b27eSBarry Smith 
10401e07b27eSBarry Smith   PetscFunctionBegin;
10411e07b27eSBarry Smith   ierr = PetscNewLog(pc,&sred);CHKERRQ(ierr);
1042*48a10b22SPatrick Sanan   sred->subcommtype    = PETSC_SUBCOMM_INTERLACED;
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);
1063*48a10b22SPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetSubcommType_C",PCTelescopeGetSubcommType_Telescope);CHKERRQ(ierr);
1064*48a10b22SPatrick 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