xref: /petsc/src/ksp/pc/impls/telescope/telescope.c (revision d8b9d5b72b0f0a8eb2091f68bce2d33abd724c2f)
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 
81e07b27eSBarry Smith /*
91e07b27eSBarry Smith  PCTelescopeSetUp_default()
101e07b27eSBarry Smith  PCTelescopeMatCreate_default()
111e07b27eSBarry Smith 
121e07b27eSBarry Smith  default
131e07b27eSBarry Smith 
141e07b27eSBarry Smith  // scatter in
151e07b27eSBarry Smith  x(comm) -> xtmp(comm)
161e07b27eSBarry Smith 
171e07b27eSBarry Smith  xred(subcomm) <- xtmp
181e07b27eSBarry Smith  yred(subcomm)
191e07b27eSBarry Smith 
201e07b27eSBarry Smith  yred(subcomm) --> xtmp
211e07b27eSBarry Smith 
221e07b27eSBarry Smith  // scatter out
231e07b27eSBarry Smith  xtmp(comm) -> y(comm)
241e07b27eSBarry Smith */
251e07b27eSBarry Smith 
261e07b27eSBarry Smith PetscBool isActiveRank(PetscSubcomm scomm)
271e07b27eSBarry Smith {
281e07b27eSBarry Smith   if (scomm->color == 0) { return PETSC_TRUE; }
291e07b27eSBarry Smith   else { return PETSC_FALSE; }
301e07b27eSBarry Smith }
311e07b27eSBarry Smith 
321e07b27eSBarry Smith #undef __FUNCT__
331e07b27eSBarry Smith #define __FUNCT__ "private_PCTelescopeGetSubDM"
341e07b27eSBarry Smith DM private_PCTelescopeGetSubDM(PC_Telescope sred)
351e07b27eSBarry Smith {
36c6a0d831SBarry Smith   DM subdm = NULL;
371e07b27eSBarry Smith 
381e07b27eSBarry Smith   if (!isActiveRank(sred->psubcomm)) { subdm = NULL; }
391e07b27eSBarry Smith   else {
401e07b27eSBarry Smith     switch (sred->sr_type) {
411e07b27eSBarry Smith     case TELESCOPE_DEFAULT: subdm = NULL;
421e07b27eSBarry Smith       break;
431e07b27eSBarry Smith     case TELESCOPE_DMDA:    subdm = ((PC_Telescope_DMDACtx*)sred->dm_ctx)->dmrepart;
441e07b27eSBarry Smith       break;
451e07b27eSBarry Smith     case TELESCOPE_DMPLEX:  subdm = NULL;
461e07b27eSBarry Smith       break;
471e07b27eSBarry Smith     }
481e07b27eSBarry Smith   }
491e07b27eSBarry Smith   return(subdm);
501e07b27eSBarry Smith }
511e07b27eSBarry Smith 
521e07b27eSBarry Smith #undef __FUNCT__
531e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeSetUp_default"
541e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_default(PC pc,PC_Telescope sred)
551e07b27eSBarry Smith {
561e07b27eSBarry Smith   PetscErrorCode ierr;
571e07b27eSBarry Smith   PetscInt       m,M,bs,st,ed;
581e07b27eSBarry Smith   Vec            x,xred,yred,xtmp;
591e07b27eSBarry Smith   Mat            B;
601e07b27eSBarry Smith   MPI_Comm       comm,subcomm;
611e07b27eSBarry Smith   VecScatter     scatter;
621e07b27eSBarry Smith   IS             isin;
631e07b27eSBarry Smith 
641e07b27eSBarry Smith   PetscFunctionBegin;
651e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: setup (default)\n");CHKERRQ(ierr);
661e07b27eSBarry Smith   comm = PetscSubcommParent(sred->psubcomm);
671e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
681e07b27eSBarry Smith 
691e07b27eSBarry Smith   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
701e07b27eSBarry Smith   ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr);
711e07b27eSBarry Smith   ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr);
721e07b27eSBarry Smith   ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr);
731e07b27eSBarry Smith 
741e07b27eSBarry Smith   xred = NULL;
753ac26c5eSBarry Smith   m    = 0;
761e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
771e07b27eSBarry Smith     ierr = VecCreate(subcomm,&xred);CHKERRQ(ierr);
781e07b27eSBarry Smith     ierr = VecSetSizes(xred,PETSC_DECIDE,M);CHKERRQ(ierr);
791e07b27eSBarry Smith     ierr = VecSetBlockSize(xred,bs);CHKERRQ(ierr);
801e07b27eSBarry Smith     ierr = VecSetFromOptions(xred);CHKERRQ(ierr);
81ca43db0aSBarry Smith     ierr = VecGetLocalSize(xred,&m);CHKERRQ(ierr);
821e07b27eSBarry Smith   }
831e07b27eSBarry Smith 
841e07b27eSBarry Smith   yred = NULL;
851e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
861e07b27eSBarry Smith     ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr);
871e07b27eSBarry Smith   }
881e07b27eSBarry Smith 
891e07b27eSBarry Smith   ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr);
901e07b27eSBarry Smith   ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr);
911e07b27eSBarry Smith   ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr);
921e07b27eSBarry Smith   ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr);
931e07b27eSBarry Smith 
941e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
951e07b27eSBarry Smith     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
961e07b27eSBarry Smith     ierr = ISCreateStride(comm,(ed-st),st,1,&isin);CHKERRQ(ierr);
971e07b27eSBarry Smith   } else {
981e07b27eSBarry Smith     ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr);
993ac26c5eSBarry Smith     ierr = ISCreateStride(comm,0,st,1,&isin);CHKERRQ(ierr);
1001e07b27eSBarry Smith   }
1011e07b27eSBarry Smith   ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr);
1021e07b27eSBarry Smith 
1031e07b27eSBarry Smith   ierr = VecScatterCreate(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr);
1041e07b27eSBarry Smith 
1051e07b27eSBarry Smith   sred->isin    = isin;
1061e07b27eSBarry Smith   sred->scatter = scatter;
1071e07b27eSBarry Smith   sred->xred    = xred;
1081e07b27eSBarry Smith   sred->yred    = yred;
1091e07b27eSBarry Smith   sred->xtmp    = xtmp;
1101e07b27eSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
1111e07b27eSBarry Smith   PetscFunctionReturn(0);
1121e07b27eSBarry Smith }
1131e07b27eSBarry Smith 
1141e07b27eSBarry Smith #undef __FUNCT__
1151e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeMatCreate_default"
1161e07b27eSBarry Smith PetscErrorCode PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A)
1171e07b27eSBarry Smith {
1181e07b27eSBarry Smith   PetscErrorCode ierr;
1191e07b27eSBarry Smith   MPI_Comm       comm,subcomm;
1201e07b27eSBarry Smith   Mat            Bred,B;
1211e07b27eSBarry Smith   PetscInt       nr,nc;
1221e07b27eSBarry Smith   IS             isrow,iscol;
1231e07b27eSBarry Smith   Mat            Blocal,*_Blocal;
1241e07b27eSBarry Smith 
1251e07b27eSBarry Smith   PetscFunctionBegin;
1261e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (default)\n");CHKERRQ(ierr);
1271e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1281e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
1291e07b27eSBarry Smith   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
1301e07b27eSBarry Smith   ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr);
1311e07b27eSBarry Smith   isrow = sred->isin;
1321e07b27eSBarry Smith   ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr);
1331e07b27eSBarry Smith   ierr = MatGetSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr);
1341e07b27eSBarry Smith   Blocal = *_Blocal;
1351e07b27eSBarry Smith   ierr = PetscFree(_Blocal);CHKERRQ(ierr);
1361e07b27eSBarry Smith   Bred = NULL;
1371e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
1381e07b27eSBarry Smith     PetscInt mm;
1391e07b27eSBarry Smith 
1401e07b27eSBarry Smith     if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; }
1411e07b27eSBarry Smith 
1421e07b27eSBarry Smith     ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr);
1431e07b27eSBarry Smith     ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr);
1441e07b27eSBarry Smith   }
1451e07b27eSBarry Smith   *A = Bred;
1461e07b27eSBarry Smith   ierr = ISDestroy(&iscol);CHKERRQ(ierr);
1471e07b27eSBarry Smith   ierr = MatDestroy(&Blocal);CHKERRQ(ierr);
1481e07b27eSBarry Smith   PetscFunctionReturn(0);
1491e07b27eSBarry Smith }
1501e07b27eSBarry Smith 
1511e07b27eSBarry Smith #undef __FUNCT__
1521e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeMatNullSpaceCreate_default"
153*d8b9d5b7SPatrick Sanan static PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat,PetscBool near)
1541e07b27eSBarry Smith {
1551e07b27eSBarry Smith   PetscErrorCode   ierr;
1561e07b27eSBarry Smith   MatNullSpace     nullspace,sub_nullspace;
1571e07b27eSBarry Smith   Mat              A,B;
1581e07b27eSBarry Smith   PetscBool        has_const;
159a947c41eSDave May   PetscInt         i,k,n = 0;
1601e07b27eSBarry Smith   const Vec        *vecs;
161c41e779fSDave May   Vec              *sub_vecs = NULL;
1621e07b27eSBarry Smith   MPI_Comm         subcomm;
1631e07b27eSBarry Smith 
1641e07b27eSBarry Smith   PetscFunctionBegin;
1651e07b27eSBarry Smith   ierr = PCGetOperators(pc,&A,&B);CHKERRQ(ierr);
166*d8b9d5b7SPatrick Sanan 
167*d8b9d5b7SPatrick Sanan   /* If the "near" flag was provided, instead operate on the near nullspace */
168*d8b9d5b7SPatrick Sanan   if (near) {
169*d8b9d5b7SPatrick Sanan     ierr = MatGetNearNullSpace(B,&nullspace);CHKERRQ(ierr);
170*d8b9d5b7SPatrick Sanan   } else {
1711e07b27eSBarry Smith     ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr);
172*d8b9d5b7SPatrick Sanan   }
173a947c41eSDave May   if (!nullspace) PetscFunctionReturn(0);
1741e07b27eSBarry Smith 
175*d8b9d5b7SPatrick Sanan   if (near) {
176*d8b9d5b7SPatrick Sanan     ierr = PetscInfo(pc,"PCTelescope: generating near nullspace (default)\n");CHKERRQ(ierr);
177*d8b9d5b7SPatrick Sanan   } else {
1781e07b27eSBarry Smith     ierr = PetscInfo(pc,"PCTelescope: generating nullspace (default)\n");CHKERRQ(ierr);
179*d8b9d5b7SPatrick Sanan   }
1801e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
1811e07b27eSBarry Smith   ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr);
1821e07b27eSBarry Smith 
1831e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
184e3acf2f7SBarry Smith     if (n) {
185e3acf2f7SBarry Smith       ierr = VecDuplicateVecs(sred->xred,n,&sub_vecs);CHKERRQ(ierr);
1861e07b27eSBarry Smith     }
1871e07b27eSBarry Smith   }
1881e07b27eSBarry Smith 
1891e07b27eSBarry Smith   /* copy entries */
1901e07b27eSBarry Smith   for (k=0; k<n; k++) {
1911e07b27eSBarry Smith     const PetscScalar *x_array;
1921e07b27eSBarry Smith     PetscScalar       *LA_sub_vec;
19313c30530SDave May     PetscInt          st,ed;
1941e07b27eSBarry Smith 
1951e07b27eSBarry Smith     /* pull in vector x->xtmp */
1961e07b27eSBarry Smith     ierr = VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1971e07b27eSBarry Smith     ierr = VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
19847856c66SBarry Smith     if (sub_vecs) {
199a04a6428SPatrick Sanan       /* copy vector entries into xred */
2001e07b27eSBarry Smith       ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
201ea2b237eSDave May       if (sub_vecs[k]) {
2021e07b27eSBarry Smith         ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr);
2031e07b27eSBarry Smith         ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
2041e07b27eSBarry Smith         for (i=0; i<ed-st; i++) {
2051e07b27eSBarry Smith           LA_sub_vec[i] = x_array[i];
2061e07b27eSBarry Smith         }
2071e07b27eSBarry Smith         ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
2081e07b27eSBarry Smith       }
2091e07b27eSBarry Smith       ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
2101e07b27eSBarry Smith     }
21147856c66SBarry Smith   }
2121e07b27eSBarry Smith 
2131e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
214*d8b9d5b7SPatrick Sanan     /* create new (near) nullspace for redundant object */
2151e07b27eSBarry Smith     ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,&sub_nullspace);CHKERRQ(ierr);
216*d8b9d5b7SPatrick Sanan     if (nullspace->remove) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope");
217*d8b9d5b7SPatrick 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");
218120bdd93SDave May 
219*d8b9d5b7SPatrick Sanan     /* attach redundant (near) nullspace to Bred */
220*d8b9d5b7SPatrick Sanan     if (near) {
221*d8b9d5b7SPatrick Sanan       ierr = MatSetNearNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr);
222*d8b9d5b7SPatrick Sanan     } else {
2231e07b27eSBarry Smith       ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr);
224*d8b9d5b7SPatrick Sanan     }
225e3acf2f7SBarry Smith     ierr = VecDestroyVecs(n,&sub_vecs);CHKERRQ(ierr);
2261e07b27eSBarry Smith   }
2271e07b27eSBarry Smith   PetscFunctionReturn(0);
2281e07b27eSBarry Smith }
2291e07b27eSBarry Smith 
2301e07b27eSBarry Smith #undef __FUNCT__
2311e07b27eSBarry Smith #define __FUNCT__ "PCView_Telescope"
2321e07b27eSBarry Smith static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer)
2331e07b27eSBarry Smith {
2341e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
2351e07b27eSBarry Smith   PetscErrorCode   ierr;
2361e07b27eSBarry Smith   PetscBool        iascii,isstring;
2371e07b27eSBarry Smith   PetscViewer      subviewer;
2381e07b27eSBarry Smith 
2391e07b27eSBarry Smith   PetscFunctionBegin;
2401e07b27eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
2411e07b27eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
2421e07b27eSBarry Smith   if (iascii) {
2431e07b27eSBarry Smith     if (!sred->psubcomm) {
2441e07b27eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: preconditioner not yet setup\n");CHKERRQ(ierr);
2451e07b27eSBarry Smith     } else {
2461e07b27eSBarry Smith       MPI_Comm    comm,subcomm;
2471e07b27eSBarry Smith       PetscMPIInt comm_size,subcomm_size;
2481e07b27eSBarry Smith       DM          dm,subdm;
2491e07b27eSBarry Smith 
2501e07b27eSBarry Smith       ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
2511e07b27eSBarry Smith       subdm = private_PCTelescopeGetSubDM(sred);
2521e07b27eSBarry Smith       comm = PetscSubcommParent(sred->psubcomm);
2531e07b27eSBarry Smith       subcomm = PetscSubcommChild(sred->psubcomm);
2541e07b27eSBarry Smith       ierr = MPI_Comm_size(comm,&comm_size);CHKERRQ(ierr);
2551e07b27eSBarry Smith       ierr = MPI_Comm_size(subcomm,&subcomm_size);CHKERRQ(ierr);
2561e07b27eSBarry Smith 
2571e07b27eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: parent comm size reduction factor = %D\n",sred->redfactor);CHKERRQ(ierr);
2581e07b27eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: comm_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size);CHKERRQ(ierr);
25948a10b22SPatrick Sanan       switch (sred->subcommtype) {
26048a10b22SPatrick Sanan         case PETSC_SUBCOMM_INTERLACED :
26148a10b22SPatrick Sanan           ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: subcomm type: interlaced\n",sred->subcommtype);CHKERRQ(ierr);
26248a10b22SPatrick Sanan           break;
26348a10b22SPatrick Sanan         case PETSC_SUBCOMM_CONTIGUOUS :
26448a10b22SPatrick Sanan           ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: subcomm type: contiguous\n",sred->subcommtype);CHKERRQ(ierr);
26548a10b22SPatrick Sanan           break;
26648a10b22SPatrick Sanan         default :
26748a10b22SPatrick Sanan           SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"General subcomm type not supported by PCTelescope");
26848a10b22SPatrick Sanan       }
2691e07b27eSBarry Smith       ierr = PetscViewerGetSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
2701e07b27eSBarry Smith       if (isActiveRank(sred->psubcomm)) {
2711e07b27eSBarry Smith         ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr);
2721e07b27eSBarry Smith 
2731e07b27eSBarry Smith         if (dm && sred->ignore_dm) {
2741e07b27eSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: ignoring DM\n");CHKERRQ(ierr);
2751e07b27eSBarry Smith         }
2767c5279cbSDave May         if (sred->ignore_kspcomputeoperators) {
2777c5279cbSDave May           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: ignoring KSPComputeOperators\n");CHKERRQ(ierr);
2787c5279cbSDave May         }
2791e07b27eSBarry Smith         switch (sred->sr_type) {
2801e07b27eSBarry Smith         case TELESCOPE_DEFAULT:
2811e07b27eSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: using default setup\n");CHKERRQ(ierr);
2821e07b27eSBarry Smith           break;
2831e07b27eSBarry Smith         case TELESCOPE_DMDA:
2841e07b27eSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: DMDA detected\n");CHKERRQ(ierr);
2851e07b27eSBarry Smith           ierr = DMView_DMDAShort(subdm,subviewer);CHKERRQ(ierr);
2861e07b27eSBarry Smith           break;
2871e07b27eSBarry Smith         case TELESCOPE_DMPLEX:
2881e07b27eSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: DMPLEX detected\n");CHKERRQ(ierr);
2891e07b27eSBarry Smith           break;
2901e07b27eSBarry Smith         }
2911e07b27eSBarry Smith 
2921e07b27eSBarry Smith         ierr = KSPView(sred->ksp,subviewer);CHKERRQ(ierr);
2931e07b27eSBarry Smith         ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr);
2941e07b27eSBarry Smith       }
2951e07b27eSBarry Smith       ierr = PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
2961e07b27eSBarry Smith     }
2971e07b27eSBarry Smith   }
2981e07b27eSBarry Smith   PetscFunctionReturn(0);
2991e07b27eSBarry Smith }
3001e07b27eSBarry Smith 
3011e07b27eSBarry Smith #undef __FUNCT__
3021e07b27eSBarry Smith #define __FUNCT__ "PCSetUp_Telescope"
3031e07b27eSBarry Smith static PetscErrorCode PCSetUp_Telescope(PC pc)
3041e07b27eSBarry Smith {
3051e07b27eSBarry Smith   PC_Telescope      sred = (PC_Telescope)pc->data;
3061e07b27eSBarry Smith   PetscErrorCode    ierr;
307bd49479cSSatish Balay   MPI_Comm          comm,subcomm=0;
3081e07b27eSBarry Smith   PCTelescopeType   sr_type;
3091e07b27eSBarry Smith 
3101e07b27eSBarry Smith   PetscFunctionBegin;
3111e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
3121e07b27eSBarry Smith 
3131e07b27eSBarry Smith   /* subcomm definition */
3141e07b27eSBarry Smith   if (!pc->setupcalled) {
3151e07b27eSBarry Smith     if (!sred->psubcomm) {
3161e07b27eSBarry Smith       ierr = PetscSubcommCreate(comm,&sred->psubcomm);CHKERRQ(ierr);
3171e07b27eSBarry Smith       ierr = PetscSubcommSetNumber(sred->psubcomm,sred->redfactor);CHKERRQ(ierr);
31848a10b22SPatrick Sanan       ierr = PetscSubcommSetType(sred->psubcomm,sred->subcommtype);CHKERRQ(ierr);
3191e07b27eSBarry Smith       ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr);
3201e07b27eSBarry Smith     }
3211e07b27eSBarry Smith   }
3220f6f40a7SSatish Balay   subcomm = PetscSubcommChild(sred->psubcomm);
3231e07b27eSBarry Smith 
3241e07b27eSBarry Smith   /* internal KSP */
3251e07b27eSBarry Smith   if (!pc->setupcalled) {
3261e07b27eSBarry Smith     const char *prefix;
3271e07b27eSBarry Smith 
3281e07b27eSBarry Smith     if (isActiveRank(sred->psubcomm)) {
3291e07b27eSBarry Smith       ierr = KSPCreate(subcomm,&sred->ksp);CHKERRQ(ierr);
3301e07b27eSBarry Smith       ierr = KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure);CHKERRQ(ierr);
3311e07b27eSBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
3321e07b27eSBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp);CHKERRQ(ierr);
3331e07b27eSBarry Smith       ierr = KSPSetType(sred->ksp,KSPPREONLY);CHKERRQ(ierr);
3341e07b27eSBarry Smith       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
3351e07b27eSBarry Smith       ierr = KSPSetOptionsPrefix(sred->ksp,prefix);CHKERRQ(ierr);
3361e07b27eSBarry Smith       ierr = KSPAppendOptionsPrefix(sred->ksp,"telescope_");CHKERRQ(ierr);
3371e07b27eSBarry Smith     }
3381e07b27eSBarry Smith   }
3391e07b27eSBarry Smith   /* Determine type of setup/update */
3401e07b27eSBarry Smith   if (!pc->setupcalled) {
3411e07b27eSBarry Smith     PetscBool has_dm,same;
3421e07b27eSBarry Smith     DM        dm;
3431e07b27eSBarry Smith 
3441e07b27eSBarry Smith     sr_type = TELESCOPE_DEFAULT;
3451e07b27eSBarry Smith     has_dm = PETSC_FALSE;
3461e07b27eSBarry Smith     ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
3471e07b27eSBarry Smith     if (dm) { has_dm = PETSC_TRUE; }
3481e07b27eSBarry Smith     if (has_dm) {
3491e07b27eSBarry Smith       /* check for dmda */
3501e07b27eSBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&same);CHKERRQ(ierr);
3511e07b27eSBarry Smith       if (same) {
3521e07b27eSBarry Smith         ierr = PetscInfo(pc,"PCTelescope: found DMDA\n");CHKERRQ(ierr);
3531e07b27eSBarry Smith         sr_type = TELESCOPE_DMDA;
3541e07b27eSBarry Smith       }
3551e07b27eSBarry Smith       /* check for dmplex */
3561e07b27eSBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same);CHKERRQ(ierr);
3571e07b27eSBarry Smith       if (same) {
3581e07b27eSBarry Smith         PetscInfo(pc,"PCTelescope: found DMPLEX\n");
3591e07b27eSBarry Smith         sr_type = TELESCOPE_DMPLEX;
3601e07b27eSBarry Smith       }
3611e07b27eSBarry Smith     }
3621e07b27eSBarry Smith 
3631e07b27eSBarry Smith     if (sred->ignore_dm) {
3641e07b27eSBarry Smith       PetscInfo(pc,"PCTelescope: ignore DM\n");
3651e07b27eSBarry Smith       sr_type = TELESCOPE_DEFAULT;
3661e07b27eSBarry Smith     }
3671e07b27eSBarry Smith     sred->sr_type = sr_type;
3681e07b27eSBarry Smith   } else {
3691e07b27eSBarry Smith     sr_type = sred->sr_type;
3701e07b27eSBarry Smith   }
3711e07b27eSBarry Smith 
372*d8b9d5b7SPatrick Sanan   /* set function pointers for repartition setup, matrix creation/update, matrix (near) nullspace, and reset functionality */
3731e07b27eSBarry Smith   switch (sr_type) {
3741e07b27eSBarry Smith   case TELESCOPE_DEFAULT:
3751e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
3761e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
3771e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
3781e07b27eSBarry Smith     sred->pctelescope_reset_type              = NULL;
3791e07b27eSBarry Smith     break;
3801e07b27eSBarry Smith   case TELESCOPE_DMDA:
3811e07b27eSBarry Smith     pc->ops->apply                            = PCApply_Telescope_dmda;
382f650675bSDave May     pc->ops->applyrichardson                  = PCApplyRichardson_Telescope_dmda;
3831e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_dmda;
3841e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_dmda;
3851e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda;
3861e07b27eSBarry Smith     sred->pctelescope_reset_type              = PCReset_Telescope_dmda;
3871e07b27eSBarry Smith     break;
388*d8b9d5b7SPatrick Sanan   case TELESCOPE_DMPLEX: SETERRQ(comm,PETSC_ERR_SUP,"Support for DMPLEX is currently not available");
3891e07b27eSBarry Smith     break;
390a04a6428SPatrick Sanan   default: SETERRQ(comm,PETSC_ERR_SUP,"Only support for repartitioning DMDA is provided");
3911e07b27eSBarry Smith     break;
3921e07b27eSBarry Smith   }
3931e07b27eSBarry Smith 
3941e07b27eSBarry Smith   /* setup */
3951e07b27eSBarry Smith   if (sred->pctelescope_setup_type) {
3961e07b27eSBarry Smith     ierr = sred->pctelescope_setup_type(pc,sred);CHKERRQ(ierr);
3971e07b27eSBarry Smith   }
3981e07b27eSBarry Smith   /* update */
3991e07b27eSBarry Smith   if (!pc->setupcalled) {
4001e07b27eSBarry Smith     if (sred->pctelescope_matcreate_type) {
4011e07b27eSBarry Smith       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred);CHKERRQ(ierr);
4021e07b27eSBarry Smith     }
4031e07b27eSBarry Smith     if (sred->pctelescope_matnullspacecreate_type) {
404*d8b9d5b7SPatrick Sanan       /* Propagate attached nullspace and near nullspace, if they exist */
405*d8b9d5b7SPatrick Sanan       ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred,PETSC_TRUE);CHKERRQ(ierr);
406*d8b9d5b7SPatrick Sanan       ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred,PETSC_FALSE);CHKERRQ(ierr);
4071e07b27eSBarry Smith     }
4081e07b27eSBarry Smith   } else {
4091e07b27eSBarry Smith     if (sred->pctelescope_matcreate_type) {
4101e07b27eSBarry Smith       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred);CHKERRQ(ierr);
4111e07b27eSBarry Smith     }
4121e07b27eSBarry Smith   }
4131e07b27eSBarry Smith 
4141e07b27eSBarry Smith   /* common - no construction */
4151e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
4161e07b27eSBarry Smith     ierr = KSPSetOperators(sred->ksp,sred->Bred,sred->Bred);CHKERRQ(ierr);
4171e07b27eSBarry Smith     if (pc->setfromoptionscalled && !pc->setupcalled){
4181e07b27eSBarry Smith       ierr = KSPSetFromOptions(sred->ksp);CHKERRQ(ierr);
4191e07b27eSBarry Smith     }
4201e07b27eSBarry Smith   }
4211e07b27eSBarry Smith   PetscFunctionReturn(0);
4221e07b27eSBarry Smith }
4231e07b27eSBarry Smith 
4241e07b27eSBarry Smith #undef __FUNCT__
4251e07b27eSBarry Smith #define __FUNCT__ "PCApply_Telescope"
4261e07b27eSBarry Smith static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y)
4271e07b27eSBarry Smith {
4281e07b27eSBarry Smith   PC_Telescope      sred = (PC_Telescope)pc->data;
4291e07b27eSBarry Smith   PetscErrorCode    ierr;
4301e07b27eSBarry Smith   Vec               xtmp,xred,yred;
43113c30530SDave May   PetscInt          i,st,ed;
4321e07b27eSBarry Smith   VecScatter        scatter;
4331e07b27eSBarry Smith   PetscScalar       *array;
4341e07b27eSBarry Smith   const PetscScalar *x_array;
4351e07b27eSBarry Smith 
4361e07b27eSBarry Smith   PetscFunctionBegin;
4371e07b27eSBarry Smith   xtmp    = sred->xtmp;
4381e07b27eSBarry Smith   scatter = sred->scatter;
4391e07b27eSBarry Smith   xred    = sred->xred;
4401e07b27eSBarry Smith   yred    = sred->yred;
4411e07b27eSBarry Smith 
4421e07b27eSBarry Smith   /* pull in vector x->xtmp */
4431e07b27eSBarry Smith   ierr = VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4441e07b27eSBarry Smith   ierr = VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4451e07b27eSBarry Smith 
4461e07b27eSBarry Smith   /* copy vector entires into xred */
4471e07b27eSBarry Smith   ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr);
4481e07b27eSBarry Smith   if (xred) {
4491e07b27eSBarry Smith     PetscScalar *LA_xred;
4501e07b27eSBarry Smith     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
4511e07b27eSBarry Smith     ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr);
4521e07b27eSBarry Smith     for (i=0; i<ed-st; i++) {
4531e07b27eSBarry Smith       LA_xred[i] = x_array[i];
4541e07b27eSBarry Smith     }
4551e07b27eSBarry Smith     ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr);
4561e07b27eSBarry Smith   }
4571e07b27eSBarry Smith   ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr);
4581e07b27eSBarry Smith   /* solve */
4591e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
4601e07b27eSBarry Smith     ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr);
4611e07b27eSBarry Smith   }
4621e07b27eSBarry Smith   /* return vector */
4631e07b27eSBarry Smith   ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr);
4641e07b27eSBarry Smith   if (yred) {
4651e07b27eSBarry Smith     const PetscScalar *LA_yred;
4661e07b27eSBarry Smith     ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr);
4671e07b27eSBarry Smith     ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr);
4681e07b27eSBarry Smith     for (i=0; i<ed-st; i++) {
4691e07b27eSBarry Smith       array[i] = LA_yred[i];
4701e07b27eSBarry Smith     }
4711e07b27eSBarry Smith     ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr);
4721e07b27eSBarry Smith   }
4731e07b27eSBarry Smith   ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr);
4741e07b27eSBarry Smith   ierr = VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4751e07b27eSBarry Smith   ierr = VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4761e07b27eSBarry Smith   PetscFunctionReturn(0);
4771e07b27eSBarry Smith }
4781e07b27eSBarry Smith 
4791e07b27eSBarry Smith #undef __FUNCT__
480f650675bSDave May #define __FUNCT__ "PCApplyRichardson_Telescope"
481f650675bSDave 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)
482f650675bSDave May {
483f650675bSDave May   PC_Telescope      sred = (PC_Telescope)pc->data;
484f650675bSDave May   PetscErrorCode    ierr;
485a1d91a28SDave May   Vec               xtmp,yred;
486f650675bSDave May   PetscInt          i,st,ed;
487f650675bSDave May   VecScatter        scatter;
488f650675bSDave May   const PetscScalar *x_array;
489f650675bSDave May   PetscBool         default_init_guess_value;
490f650675bSDave May 
491f650675bSDave May   PetscFunctionBegin;
492f650675bSDave May   xtmp    = sred->xtmp;
493f650675bSDave May   scatter = sred->scatter;
494f650675bSDave May   yred    = sred->yred;
495f650675bSDave May 
496f650675bSDave May   if (its > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCApplyRichardson_Telescope only supports max_it = 1");
497f650675bSDave May   *reason = (PCRichardsonConvergedReason)0;
498f650675bSDave May 
499f650675bSDave May   if (!zeroguess) {
500f650675bSDave May     ierr = PetscInfo(pc,"PCTelescope: Scattering y for non-zero initial guess\n");CHKERRQ(ierr);
501f650675bSDave May     /* pull in vector y->xtmp */
502f650675bSDave May     ierr = VecScatterBegin(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
503f650675bSDave May     ierr = VecScatterEnd(scatter,y,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
504f650675bSDave May 
505f650675bSDave May     /* copy vector entires into xred */
506f650675bSDave May     ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr);
507f650675bSDave May     if (yred) {
508f650675bSDave May       PetscScalar *LA_yred;
509f650675bSDave May       ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr);
510f650675bSDave May       ierr = VecGetArray(yred,&LA_yred);CHKERRQ(ierr);
511f650675bSDave May       for (i=0; i<ed-st; i++) {
512f650675bSDave May         LA_yred[i] = x_array[i];
513f650675bSDave May       }
514f650675bSDave May       ierr = VecRestoreArray(yred,&LA_yred);CHKERRQ(ierr);
515f650675bSDave May     }
516f650675bSDave May     ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr);
517f650675bSDave May   }
518f650675bSDave May 
519f650675bSDave May   if (isActiveRank(sred->psubcomm)) {
520f650675bSDave May     ierr = KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value);CHKERRQ(ierr);
521f650675bSDave May     if (!zeroguess) ierr = KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE);CHKERRQ(ierr);
522f650675bSDave May   }
523f650675bSDave May 
524f650675bSDave May   ierr = PCApply_Telescope(pc,x,y);CHKERRQ(ierr);
525f650675bSDave May 
526f650675bSDave May   if (isActiveRank(sred->psubcomm)) {
527f650675bSDave May     ierr = KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value);CHKERRQ(ierr);
528f650675bSDave May   }
529f650675bSDave May 
530f650675bSDave May   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
531f650675bSDave May   *outits = 1;
532f650675bSDave May   PetscFunctionReturn(0);
533f650675bSDave May }
534f650675bSDave May 
535f650675bSDave May #undef __FUNCT__
5361e07b27eSBarry Smith #define __FUNCT__ "PCReset_Telescope"
5371e07b27eSBarry Smith static PetscErrorCode PCReset_Telescope(PC pc)
5381e07b27eSBarry Smith {
5391e07b27eSBarry Smith   PC_Telescope   sred = (PC_Telescope)pc->data;
5401e07b27eSBarry Smith   PetscErrorCode ierr;
5411e07b27eSBarry Smith 
5421e07b27eSBarry Smith   ierr = ISDestroy(&sred->isin);CHKERRQ(ierr);
5431e07b27eSBarry Smith   ierr = VecScatterDestroy(&sred->scatter);CHKERRQ(ierr);
544e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->xred);CHKERRQ(ierr);
545e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->yred);CHKERRQ(ierr);
546e3acf2f7SBarry Smith   ierr = VecDestroy(&sred->xtmp);CHKERRQ(ierr);
547e3acf2f7SBarry Smith   ierr = MatDestroy(&sred->Bred);CHKERRQ(ierr);
548e3acf2f7SBarry Smith   ierr = KSPReset(sred->ksp);CHKERRQ(ierr);
5491e07b27eSBarry Smith   if (sred->pctelescope_reset_type) {
5501e07b27eSBarry Smith     ierr = sred->pctelescope_reset_type(pc);CHKERRQ(ierr);
5511e07b27eSBarry Smith   }
5521e07b27eSBarry Smith   PetscFunctionReturn(0);
5531e07b27eSBarry Smith }
5541e07b27eSBarry Smith 
5551e07b27eSBarry Smith #undef __FUNCT__
5561e07b27eSBarry Smith #define __FUNCT__ "PCDestroy_Telescope"
5571e07b27eSBarry Smith static PetscErrorCode PCDestroy_Telescope(PC pc)
5581e07b27eSBarry Smith {
5591e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
5601e07b27eSBarry Smith   PetscErrorCode   ierr;
5611e07b27eSBarry Smith 
5621e07b27eSBarry Smith   PetscFunctionBegin;
5631e07b27eSBarry Smith   ierr = PCReset_Telescope(pc);CHKERRQ(ierr);
564e3acf2f7SBarry Smith   ierr = KSPDestroy(&sred->ksp);CHKERRQ(ierr);
5651e07b27eSBarry Smith   ierr = PetscSubcommDestroy(&sred->psubcomm);CHKERRQ(ierr);
566e3acf2f7SBarry Smith   ierr = PetscFree(sred->dm_ctx);CHKERRQ(ierr);
567e3acf2f7SBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
5681e07b27eSBarry Smith   PetscFunctionReturn(0);
5691e07b27eSBarry Smith }
5701e07b27eSBarry Smith 
5711e07b27eSBarry Smith #undef __FUNCT__
5721e07b27eSBarry Smith #define __FUNCT__ "PCSetFromOptions_Telescope"
5734416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Telescope(PetscOptionItems *PetscOptionsObject,PC pc)
5741e07b27eSBarry Smith {
5751e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
5761e07b27eSBarry Smith   PetscErrorCode   ierr;
5771e07b27eSBarry Smith   MPI_Comm         comm;
5781e07b27eSBarry Smith   PetscMPIInt      size;
57948a10b22SPatrick Sanan   PetscBool        flg;
58048a10b22SPatrick Sanan   PetscSubcommType subcommtype;
5811e07b27eSBarry Smith 
5821e07b27eSBarry Smith   PetscFunctionBegin;
5831e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
5841e07b27eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
5851e07b27eSBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Telescope options");CHKERRQ(ierr);
58648a10b22SPatrick Sanan   ierr = PetscOptionsEnum("-pc_telescope_subcomm_type","Subcomm type (interlaced or contiguous)","PCTelescopeSetSubcommType",PetscSubcommTypes,(PetscEnum)sred->subcommtype,(PetscEnum*)&subcommtype,&flg);CHKERRQ(ierr);
58748a10b22SPatrick Sanan   if (flg) {
58848a10b22SPatrick Sanan     ierr = PCTelescopeSetSubcommType(pc,subcommtype);CHKERRQ(ierr);
58948a10b22SPatrick Sanan   }
5901e07b27eSBarry Smith   ierr = PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,0);CHKERRQ(ierr);
5911e07b27eSBarry Smith   if (sred->redfactor > size) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size");
5921e07b27eSBarry Smith   ierr = PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,0);CHKERRQ(ierr);
5937c5279cbSDave May   ierr = PetscOptionsBool("-pc_telescope_ignore_kspcomputeoperators","Ignore method used to compute A","PCTelescopeSetIgnoreKSPComputeOperators",sred->ignore_kspcomputeoperators,&sred->ignore_kspcomputeoperators,0);CHKERRQ(ierr);
5941e07b27eSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
5951e07b27eSBarry Smith   PetscFunctionReturn(0);
5961e07b27eSBarry Smith }
5971e07b27eSBarry Smith 
5981e07b27eSBarry Smith /* PC simplementation specific API's */
5991e07b27eSBarry Smith 
60048a10b22SPatrick Sanan #undef __FUNCT__
60148a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetKSP_Telescope"
6021e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp)
6031e07b27eSBarry Smith {
6041e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
605bd49479cSSatish Balay   PetscFunctionBegin;
6061e07b27eSBarry Smith   if (ksp) *ksp = red->ksp;
607bd49479cSSatish Balay   PetscFunctionReturn(0);
6081e07b27eSBarry Smith }
6091e07b27eSBarry Smith 
61048a10b22SPatrick Sanan #undef __FUNCT__
61148a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetSubcommType_Telescope"
61248a10b22SPatrick Sanan static PetscErrorCode PCTelescopeGetSubcommType_Telescope(PC pc,PetscSubcommType *subcommtype)
61348a10b22SPatrick Sanan {
61448a10b22SPatrick Sanan   PC_Telescope red = (PC_Telescope)pc->data;
61548a10b22SPatrick Sanan   PetscFunctionBegin;
61648a10b22SPatrick Sanan   if (subcommtype) *subcommtype = red->subcommtype;
61748a10b22SPatrick Sanan   PetscFunctionReturn(0);
61848a10b22SPatrick Sanan }
61948a10b22SPatrick Sanan 
62048a10b22SPatrick Sanan #undef __FUNCT__
62148a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeSetSubcommType_Telescope"
62248a10b22SPatrick Sanan static PetscErrorCode PCTelescopeSetSubcommType_Telescope(PC pc,PetscSubcommType subcommtype)
62348a10b22SPatrick Sanan {
62448a10b22SPatrick Sanan   PC_Telescope     red = (PC_Telescope)pc->data;
62548a10b22SPatrick Sanan 
62648a10b22SPatrick Sanan   PetscFunctionBegin;
62748a10b22SPatrick 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.");
62848a10b22SPatrick Sanan   red->subcommtype = subcommtype;
62948a10b22SPatrick Sanan   PetscFunctionReturn(0);
63048a10b22SPatrick Sanan }
63148a10b22SPatrick Sanan 
63248a10b22SPatrick Sanan #undef __FUNCT__
63348a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetReductionFactor_Telescope"
6341e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact)
6351e07b27eSBarry Smith {
6361e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
637bd49479cSSatish Balay   PetscFunctionBegin;
6381e07b27eSBarry Smith   if (fact) *fact = red->redfactor;
639bd49479cSSatish Balay   PetscFunctionReturn(0);
6401e07b27eSBarry Smith }
6411e07b27eSBarry Smith 
64248a10b22SPatrick Sanan #undef __FUNCT__
64348a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeSetReductionFactor_Telescope"
6441e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact)
6451e07b27eSBarry Smith {
6461e07b27eSBarry Smith   PC_Telescope     red = (PC_Telescope)pc->data;
6471e07b27eSBarry Smith   PetscMPIInt      size;
6481e07b27eSBarry Smith   PetscErrorCode   ierr;
6491e07b27eSBarry Smith 
650bd49479cSSatish Balay   PetscFunctionBegin;
6511e07b27eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
6521e07b27eSBarry Smith   if (fact <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact);
6531e07b27eSBarry Smith   if (fact > size) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact);
6541e07b27eSBarry Smith   red->redfactor = fact;
655bd49479cSSatish Balay   PetscFunctionReturn(0);
6561e07b27eSBarry Smith }
6571e07b27eSBarry Smith 
65848a10b22SPatrick Sanan #undef __FUNCT__
65948a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetIgnoreDM_Telescope"
6601e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v)
6611e07b27eSBarry Smith {
6621e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
663bd49479cSSatish Balay   PetscFunctionBegin;
6641e07b27eSBarry Smith   if (v) *v = red->ignore_dm;
665bd49479cSSatish Balay   PetscFunctionReturn(0);
6661e07b27eSBarry Smith }
66748a10b22SPatrick Sanan 
66848a10b22SPatrick Sanan #undef __FUNCT__
66948a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeSetIgnoreDM_Telescope"
6701e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v)
6711e07b27eSBarry Smith {
6721e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
673bd49479cSSatish Balay   PetscFunctionBegin;
6741e07b27eSBarry Smith   red->ignore_dm = v;
675bd49479cSSatish Balay   PetscFunctionReturn(0);
6761e07b27eSBarry Smith }
6771e07b27eSBarry Smith 
67848a10b22SPatrick Sanan #undef __FUNCT__
67948a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetIgnoreKSPComputeOperators_Telescope"
6800ae7c45bSDave May static PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool *v)
6810ae7c45bSDave May {
6820ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
6830ae7c45bSDave May   PetscFunctionBegin;
6840ae7c45bSDave May   if (v) *v = red->ignore_kspcomputeoperators;
6850ae7c45bSDave May   PetscFunctionReturn(0);
6860ae7c45bSDave May }
68748a10b22SPatrick Sanan 
68848a10b22SPatrick Sanan #undef __FUNCT__
68948a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeSetIgnoreKSPComputeOperators_Telescope"
6900ae7c45bSDave May static PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators_Telescope(PC pc,PetscBool v)
6910ae7c45bSDave May {
6920ae7c45bSDave May   PC_Telescope red = (PC_Telescope)pc->data;
6930ae7c45bSDave May   PetscFunctionBegin;
6940ae7c45bSDave May   red->ignore_kspcomputeoperators = v;
6950ae7c45bSDave May   PetscFunctionReturn(0);
6960ae7c45bSDave May }
6970ae7c45bSDave May 
69848a10b22SPatrick Sanan #undef __FUNCT__
69948a10b22SPatrick Sanan #define __FUNCT__ "PCTelescopeGetDM_Telescope"
7001e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm)
7011e07b27eSBarry Smith {
7021e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
703bd49479cSSatish Balay   PetscFunctionBegin;
7041e07b27eSBarry Smith   *dm = private_PCTelescopeGetSubDM(red);
705bd49479cSSatish Balay   PetscFunctionReturn(0);
7061e07b27eSBarry Smith }
7071e07b27eSBarry Smith 
7080f43ea10SBarry Smith #undef __FUNCT__
7090f43ea10SBarry Smith #define __FUNCT__ "PCTelescopeGetKSP"
7101e07b27eSBarry Smith /*@
7111e07b27eSBarry Smith  PCTelescopeGetKSP - Gets the KSP created by the telescoping PC.
7121e07b27eSBarry Smith 
7131e07b27eSBarry Smith  Not Collective
7141e07b27eSBarry Smith 
7151e07b27eSBarry Smith  Input Parameter:
7161e07b27eSBarry Smith .  pc - the preconditioner context
7171e07b27eSBarry Smith 
7181e07b27eSBarry Smith  Output Parameter:
7191e07b27eSBarry Smith .  subksp - the KSP defined the smaller set of processes
7201e07b27eSBarry Smith 
7211e07b27eSBarry Smith  Level: advanced
7221e07b27eSBarry Smith 
7231e07b27eSBarry Smith .keywords: PC, telescopting solve
7241e07b27eSBarry Smith @*/
7251e07b27eSBarry Smith PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp)
7261e07b27eSBarry Smith {
727bd49479cSSatish Balay   PetscErrorCode ierr;
728bd49479cSSatish Balay   PetscFunctionBegin;
729163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp));CHKERRQ(ierr);
730bd49479cSSatish Balay   PetscFunctionReturn(0);
7311e07b27eSBarry Smith }
7321e07b27eSBarry Smith 
7330f43ea10SBarry Smith #undef __FUNCT__
7340f43ea10SBarry Smith #define __FUNCT__ "PCTelescopeGetReductionFactor"
7351e07b27eSBarry Smith /*@
7361e07b27eSBarry Smith  PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by.
7371e07b27eSBarry Smith 
7381e07b27eSBarry Smith  Not Collective
7391e07b27eSBarry Smith 
7401e07b27eSBarry Smith  Input Parameter:
7411e07b27eSBarry Smith .  pc - the preconditioner context
7421e07b27eSBarry Smith 
7431e07b27eSBarry Smith  Output Parameter:
7441e07b27eSBarry Smith .  fact - the reduction factor
7451e07b27eSBarry Smith 
7461e07b27eSBarry Smith  Level: advanced
7471e07b27eSBarry Smith 
7481e07b27eSBarry Smith .keywords: PC, telescoping solve
7491e07b27eSBarry Smith @*/
7501e07b27eSBarry Smith PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact)
7511e07b27eSBarry Smith {
752bd49479cSSatish Balay   PetscErrorCode ierr;
753bd49479cSSatish Balay   PetscFunctionBegin;
754163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact));CHKERRQ(ierr);
755bd49479cSSatish Balay   PetscFunctionReturn(0);
7561e07b27eSBarry Smith }
7571e07b27eSBarry Smith 
7580f43ea10SBarry Smith #undef __FUNCT__
7590f43ea10SBarry Smith #define __FUNCT__ "PCTelescopeSetReductionFactor"
7601e07b27eSBarry Smith /*@
7611e07b27eSBarry Smith  PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by.
7621e07b27eSBarry Smith 
7631e07b27eSBarry Smith  Not Collective
7641e07b27eSBarry Smith 
7651e07b27eSBarry Smith  Input Parameter:
7661e07b27eSBarry Smith .  pc - the preconditioner context
7671e07b27eSBarry Smith 
7681e07b27eSBarry Smith  Output Parameter:
7691e07b27eSBarry Smith .  fact - the reduction factor
7701e07b27eSBarry Smith 
7711e07b27eSBarry Smith  Level: advanced
7721e07b27eSBarry Smith 
7731e07b27eSBarry Smith .keywords: PC, telescoping solve
7741e07b27eSBarry Smith @*/
7751e07b27eSBarry Smith PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact)
7761e07b27eSBarry Smith {
777bd49479cSSatish Balay   PetscErrorCode ierr;
778bd49479cSSatish Balay   PetscFunctionBegin;
779bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact));CHKERRQ(ierr);
780bd49479cSSatish Balay   PetscFunctionReturn(0);
7811e07b27eSBarry Smith }
7821e07b27eSBarry Smith 
7830f43ea10SBarry Smith #undef __FUNCT__
7840f43ea10SBarry Smith #define __FUNCT__ "PCTelescopeGetIgnoreDM"
7851e07b27eSBarry Smith /*@
7861e07b27eSBarry Smith  PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used.
7871e07b27eSBarry Smith 
7881e07b27eSBarry Smith  Not Collective
7891e07b27eSBarry Smith 
7901e07b27eSBarry Smith  Input Parameter:
7911e07b27eSBarry Smith .  pc - the preconditioner context
7921e07b27eSBarry Smith 
7931e07b27eSBarry Smith  Output Parameter:
7941e07b27eSBarry Smith .  v - the flag
7951e07b27eSBarry Smith 
7961e07b27eSBarry Smith  Level: advanced
7971e07b27eSBarry Smith 
7981e07b27eSBarry Smith .keywords: PC, telescoping solve
7991e07b27eSBarry Smith @*/
8001e07b27eSBarry Smith PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v)
8011e07b27eSBarry Smith {
802bd49479cSSatish Balay   PetscErrorCode ierr;
803bd49479cSSatish Balay   PetscFunctionBegin;
804163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr);
805bd49479cSSatish Balay   PetscFunctionReturn(0);
8061e07b27eSBarry Smith }
8071e07b27eSBarry Smith 
8080f43ea10SBarry Smith #undef __FUNCT__
8090f43ea10SBarry Smith #define __FUNCT__ "PCTelescopeSetIgnoreDM"
8101e07b27eSBarry Smith /*@
8111e07b27eSBarry Smith  PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC.
8121e07b27eSBarry Smith 
8131e07b27eSBarry Smith  Not Collective
8141e07b27eSBarry Smith 
8151e07b27eSBarry Smith  Input Parameter:
8161e07b27eSBarry Smith .  pc - the preconditioner context
8171e07b27eSBarry Smith 
8181e07b27eSBarry Smith  Output Parameter:
8191e07b27eSBarry Smith .  v - Use PETSC_TRUE to ignore any DM
8201e07b27eSBarry Smith 
8211e07b27eSBarry Smith  Level: advanced
8221e07b27eSBarry Smith 
8231e07b27eSBarry Smith .keywords: PC, telescoping solve
8241e07b27eSBarry Smith @*/
825bfd6bcc6SSatish Balay PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v)
8261e07b27eSBarry Smith {
827bd49479cSSatish Balay   PetscErrorCode ierr;
828bd49479cSSatish Balay   PetscFunctionBegin;
829bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr);
830bd49479cSSatish Balay   PetscFunctionReturn(0);
8311e07b27eSBarry Smith }
8321e07b27eSBarry Smith 
8330f43ea10SBarry Smith #undef __FUNCT__
8340f43ea10SBarry Smith #define __FUNCT__ "PCTelescopeGetIgnoreKSPComputeOperators"
8351e07b27eSBarry Smith /*@
8360ae7c45bSDave May  PCTelescopeGetIgnoreKSPComputeOperators - Get the flag indicating if KSPComputeOperators will be used.
8370ae7c45bSDave May 
8380ae7c45bSDave May  Not Collective
8390ae7c45bSDave May 
8400ae7c45bSDave May  Input Parameter:
8410ae7c45bSDave May .  pc - the preconditioner context
8420ae7c45bSDave May 
8430ae7c45bSDave May  Output Parameter:
8440ae7c45bSDave May .  v - the flag
8450ae7c45bSDave May 
8460ae7c45bSDave May  Level: advanced
8470ae7c45bSDave May 
8480ae7c45bSDave May .keywords: PC, telescoping solve
8490ae7c45bSDave May @*/
8500ae7c45bSDave May PetscErrorCode PCTelescopeGetIgnoreKSPComputeOperators(PC pc,PetscBool *v)
8510ae7c45bSDave May {
8520ae7c45bSDave May   PetscErrorCode ierr;
8530ae7c45bSDave May   PetscFunctionBegin;
854163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr);
8550ae7c45bSDave May   PetscFunctionReturn(0);
8560ae7c45bSDave May }
8570ae7c45bSDave May 
8580f43ea10SBarry Smith #undef __FUNCT__
8590f43ea10SBarry Smith #define __FUNCT__ "PCTelescopeSetIgnoreKSPComputeOperators"
8600ae7c45bSDave May /*@
8610ae7c45bSDave May  PCTelescopeSetIgnoreKSPComputeOperators - Set a flag to ignore KSPComputeOperators.
8620ae7c45bSDave May 
8630ae7c45bSDave May  Not Collective
8640ae7c45bSDave May 
8650ae7c45bSDave May  Input Parameter:
8660ae7c45bSDave May .  pc - the preconditioner context
8670ae7c45bSDave May 
8680ae7c45bSDave May  Output Parameter:
869a954d8f4SDave May .  v - Use PETSC_TRUE to ignore the method (if defined) set via KSPSetComputeOperators on pc
8700ae7c45bSDave May 
8710ae7c45bSDave May  Level: advanced
8720ae7c45bSDave May 
8730ae7c45bSDave May .keywords: PC, telescoping solve
8740ae7c45bSDave May @*/
8750ae7c45bSDave May PetscErrorCode PCTelescopeSetIgnoreKSPComputeOperators(PC pc,PetscBool v)
8760ae7c45bSDave May {
8770ae7c45bSDave May   PetscErrorCode ierr;
8780ae7c45bSDave May   PetscFunctionBegin;
8790ae7c45bSDave May   ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr);
8800ae7c45bSDave May   PetscFunctionReturn(0);
8810ae7c45bSDave May }
8820ae7c45bSDave May 
8830f43ea10SBarry Smith #undef __FUNCT__
8840f43ea10SBarry Smith #define __FUNCT__ "PCTelescopeGetDM"
8850ae7c45bSDave May /*@
8861e07b27eSBarry Smith  PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP.
8871e07b27eSBarry Smith 
8881e07b27eSBarry Smith  Not Collective
8891e07b27eSBarry Smith 
8901e07b27eSBarry Smith  Input Parameter:
8911e07b27eSBarry Smith .  pc - the preconditioner context
8921e07b27eSBarry Smith 
8931e07b27eSBarry Smith  Output Parameter:
8941e07b27eSBarry Smith .  subdm - The re-partitioned DM
8951e07b27eSBarry Smith 
8961e07b27eSBarry Smith  Level: advanced
8971e07b27eSBarry Smith 
8981e07b27eSBarry Smith .keywords: PC, telescoping solve
8991e07b27eSBarry Smith @*/
9001e07b27eSBarry Smith PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm)
9011e07b27eSBarry Smith {
902bd49479cSSatish Balay   PetscErrorCode ierr;
903bd49479cSSatish Balay   PetscFunctionBegin;
904163d334eSBarry Smith   ierr = PetscUseMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm));CHKERRQ(ierr);
905bd49479cSSatish Balay   PetscFunctionReturn(0);
9061e07b27eSBarry Smith }
9071e07b27eSBarry Smith 
9080f43ea10SBarry Smith #undef __FUNCT__
9090f43ea10SBarry Smith #define __FUNCT__ "PCTelescopeSetSubcommType"
91048a10b22SPatrick Sanan /*@
91148a10b22SPatrick Sanan  PCTelescopeSetSubcommType - set subcommunicator type (interlaced or contiguous)
91248a10b22SPatrick Sanan 
91348a10b22SPatrick Sanan  Logically Collective
91448a10b22SPatrick Sanan 
91548a10b22SPatrick Sanan  Input Parameter:
9161dae98e4SBarry Smith +  pc - the preconditioner context
9171dae98e4SBarry Smith -  subcommtype - the subcommunicator type (see PetscSubcommType)
91848a10b22SPatrick Sanan 
91948a10b22SPatrick Sanan  Level: advanced
92048a10b22SPatrick Sanan 
92148a10b22SPatrick Sanan .keywords: PC, telescoping solve
92248a10b22SPatrick Sanan 
92348a10b22SPatrick Sanan .seealso: PetscSubcommType, PetscSubcomm, PCTELESCOPE
92448a10b22SPatrick Sanan @*/
92548a10b22SPatrick Sanan PetscErrorCode PCTelescopeSetSubcommType(PC pc, PetscSubcommType subcommtype)
92648a10b22SPatrick Sanan {
92748a10b22SPatrick Sanan   PetscErrorCode ierr;
92848a10b22SPatrick Sanan   PetscFunctionBegin;
92948a10b22SPatrick Sanan   ierr = PetscTryMethod(pc,"PCTelescopeSetSubcommType_C",(PC,PetscSubcommType),(pc,subcommtype));CHKERRQ(ierr);
93048a10b22SPatrick Sanan   PetscFunctionReturn(0);
93148a10b22SPatrick Sanan }
93248a10b22SPatrick Sanan 
9330f43ea10SBarry Smith #undef __FUNCT__
9341dae98e4SBarry Smith #define __FUNCT__ "PCTelescopeGetSubcommType"
93548a10b22SPatrick Sanan /*@
93648a10b22SPatrick Sanan  PCTelescopeGetSubcommType - Get the subcommunicator type (interlaced or contiguous)
93748a10b22SPatrick Sanan 
93848a10b22SPatrick Sanan  Not Collective
93948a10b22SPatrick Sanan 
94048a10b22SPatrick Sanan  Input Parameter:
94148a10b22SPatrick Sanan .  pc - the preconditioner context
94248a10b22SPatrick Sanan 
94348a10b22SPatrick Sanan  Output Parameter:
94448a10b22SPatrick Sanan .  subcommtype - the subcommunicator type (see PetscSubcommType)
94548a10b22SPatrick Sanan 
94648a10b22SPatrick Sanan  Level: advanced
94748a10b22SPatrick Sanan 
94848a10b22SPatrick Sanan .keywords: PC, telescoping solve
94948a10b22SPatrick Sanan 
9501dae98e4SBarry Smith .seealso: PetscSubcomm, PetscSubcommType, PCTELESCOPE
95148a10b22SPatrick Sanan @*/
9521dae98e4SBarry Smith PetscErrorCode PCTelescopeGetSubcommType(PC pc, PetscSubcommType *subcommtype)
95348a10b22SPatrick Sanan {
95448a10b22SPatrick Sanan   PetscErrorCode ierr;
95548a10b22SPatrick Sanan   PetscFunctionBegin;
95648a10b22SPatrick Sanan   ierr = PetscUseMethod(pc,"PCTelescopeGetSubcommType_C",(PC,PetscSubcommType*),(pc,subcommtype));CHKERRQ(ierr);
95748a10b22SPatrick Sanan   PetscFunctionReturn(0);
95848a10b22SPatrick Sanan }
95948a10b22SPatrick Sanan 
9601e07b27eSBarry Smith /* -------------------------------------------------------------------------------------*/
9611e07b27eSBarry Smith /*MC
9621e07b27eSBarry Smith    PCTELESCOPE - Runs a KSP solver on a sub-group of processors. MPI processes not in the sub-communicator are idle during the solve.
9631e07b27eSBarry Smith 
9641e07b27eSBarry Smith    Options Database:
965a04a6428SPatrick 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.
966a04a6428SPatrick Sanan -  -pc_telescope_ignore_dm  - flag to indicate whether an attached DM should be ignored
967a04a6428SPatrick Sanan -  -pc_telescope_subcomm_type <interlaced,contiguous> - how to define the reduced communicator. see PetscSubcomm for more.
9681e07b27eSBarry Smith 
9691e07b27eSBarry Smith    Level: advanced
9701e07b27eSBarry Smith 
9711e07b27eSBarry Smith    Notes:
9726fc41876SBarry Smith    The preconditioner is deemed telescopic as it only calls KSPSolve() on a single
9736fc41876SBarry Smith    sub-communicator in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators.
974a04a6428SPatrick Sanan    This means there will be MPI processes within c which will be idle during the application of this preconditioner.
9756fc41876SBarry Smith 
9761e07b27eSBarry Smith    The default KSP is PREONLY. If a DM is attached to the PC, it is re-partitioned on the sub-communicator.
9771e07b27eSBarry 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.
9781e07b27eSBarry Smith    Currently only support for re-partitioning a DMDA is provided.
9791e07b27eSBarry Smith    Any nullspace attached to the original Bmat operator are extracted, re-partitioned and set on the repartitioned Bmat operator.
9801e07b27eSBarry Smith    KSPSetComputeOperators() is not propagated to the sub KSP.
9811e07b27eSBarry Smith    Currently there is no support for the flag -pc_use_amat
9821e07b27eSBarry Smith 
9836fc41876SBarry Smith    Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation
9846fc41876SBarry Smith    creates a child sub-communicator (c') containing less MPI processes than the original parent preconditioner (PC).
9856fc41876SBarry Smith 
9866fc41876SBarry Smith   Developer Notes:
9876fc41876SBarry Smith    During PCSetup, the B operator is scattered onto c'.
9886fc41876SBarry Smith    Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c').
9896fc41876SBarry Smith    Then KSPSolve() is executed on the c' communicator.
9906fc41876SBarry Smith 
9916fc41876SBarry Smith    The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED
992a04a6428SPatrick 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.
9936fc41876SBarry Smith 
9946fc41876SBarry Smith    The telescoping preconditioner is aware of nullspaces which are attached to the only B operator.
9956fc41876SBarry Smith    In case where B has a n nullspace attached, these nullspaces vectors are extract from B and mapped into
9966fc41876SBarry Smith    a new nullspace (defined on the sub-communicator) which is attached to B' (the B operator which was scattered to c')
9976fc41876SBarry Smith 
9986fc41876SBarry Smith    The telescoping preconditioner is aware of an attached DM. In the event that the DM is of type DMDA (2D or 3D -
9996fc41876SBarry 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
10006fc41876SBarry Smith    is attached the sub KSPSolve(). The design of telescope is such that it should be possible to extend support
10016fc41876SBarry Smith    for re-partitioning other DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs.
10026fc41876SBarry Smith 
10036fc41876SBarry Smith    By default, B' is defined by simply fusing rows from different MPI processes
10046fc41876SBarry Smith 
10056fc41876SBarry Smith    When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permuting of B
10066fc41876SBarry Smith    into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatGetSubMatrices(), (iii) fusing the
10076fc41876SBarry Smith    locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat()
10086fc41876SBarry Smith 
10096fc41876SBarry Smith    Limitations/improvements
10106fc41876SBarry Smith    VecPlaceArray could be used within PCApply() to improve efficiency and reduce memory usage.
10116fc41876SBarry Smith 
10126fc41876SBarry Smith    The symmetric permutation used when a DMDA is encountered is performed via explicitly assmbleming a permutation matrix P,
10136fc41876SBarry 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
10146fc41876SBarry Smith    VecPermute() does not supported for the use case required here. By computing P, I can permute both the operator and RHS in a
10156fc41876SBarry Smith    consistent manner.
10166fc41876SBarry Smith 
10176fc41876SBarry Smith    Mapping of vectors is performed this way
10186fc41876SBarry 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.
10196fc41876SBarry Smith    Using the interlaced creation routine, the ranks in c with color = 0, will be rank 0 and 2.
10206fc41876SBarry Smith    We perform the scatter to the sub-comm in the following way,
10216fc41876SBarry Smith    [1] Given a vector x defined on comm c
10226fc41876SBarry Smith 
10236fc41876SBarry Smith    rank(c) : _________ 0 ______  ________ 1 _______  ________ 2 _____________ ___________ 3 __________
10246fc41876SBarry 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]
10256fc41876SBarry Smith 
10266fc41876SBarry Smith    scatter to xtmp defined also on comm c so that we have the following values
10276fc41876SBarry Smith 
10286fc41876SBarry Smith    rank(c) : ___________________ 0 ________________  _1_  ______________________ 2 _______________________  __3_
10296fc41876SBarry 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] [  ]
10306fc41876SBarry Smith 
10316fc41876SBarry Smith    The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') have no values
10326fc41876SBarry Smith 
10336fc41876SBarry Smith 
10346fc41876SBarry 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'.
10356fc41876SBarry Smith    Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0.
10366fc41876SBarry Smith 
10376fc41876SBarry Smith     rank(c') : ___________________ 0 _______________  ______________________ 1 _____________________
10386fc41876SBarry 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]
10396fc41876SBarry Smith 
10401e07b27eSBarry Smith 
10411e07b27eSBarry Smith   Contributed by Dave May
10421e07b27eSBarry Smith 
10436fc41876SBarry Smith .seealso:  PCTelescopeGetKSP(), PCTelescopeGetDM(), PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(), PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT
10441e07b27eSBarry Smith M*/
10451e07b27eSBarry Smith #undef __FUNCT__
10461e07b27eSBarry Smith #define __FUNCT__ "PCCreate_Telescope"
10471e07b27eSBarry Smith PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc)
10481e07b27eSBarry Smith {
10491e07b27eSBarry Smith   PetscErrorCode       ierr;
10501e07b27eSBarry Smith   struct _PC_Telescope *sred;
10511e07b27eSBarry Smith 
10521e07b27eSBarry Smith   PetscFunctionBegin;
10531e07b27eSBarry Smith   ierr = PetscNewLog(pc,&sred);CHKERRQ(ierr);
105448a10b22SPatrick Sanan   sred->subcommtype    = PETSC_SUBCOMM_INTERLACED;
10551e07b27eSBarry Smith   sred->redfactor      = 1;
10561e07b27eSBarry Smith   sred->ignore_dm      = PETSC_FALSE;
10577c5279cbSDave May   sred->ignore_kspcomputeoperators = PETSC_FALSE;
10581e07b27eSBarry Smith   pc->data             = (void*)sred;
10591e07b27eSBarry Smith 
10601e07b27eSBarry Smith   pc->ops->apply           = PCApply_Telescope;
10611e07b27eSBarry Smith   pc->ops->applytranspose  = NULL;
1062f650675bSDave May   pc->ops->applyrichardson = PCApplyRichardson_Telescope;
10631e07b27eSBarry Smith   pc->ops->setup           = PCSetUp_Telescope;
10641e07b27eSBarry Smith   pc->ops->destroy         = PCDestroy_Telescope;
10651e07b27eSBarry Smith   pc->ops->reset           = PCReset_Telescope;
10661e07b27eSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Telescope;
10671e07b27eSBarry Smith   pc->ops->view            = PCView_Telescope;
10681e07b27eSBarry Smith 
10691e07b27eSBarry Smith   sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
10701e07b27eSBarry Smith   sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
10711e07b27eSBarry Smith   sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
10721e07b27eSBarry Smith   sred->pctelescope_reset_type              = NULL;
10731e07b27eSBarry Smith 
10741e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetKSP_C",PCTelescopeGetKSP_Telescope);CHKERRQ(ierr);
107548a10b22SPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetSubcommType_C",PCTelescopeGetSubcommType_Telescope);CHKERRQ(ierr);
107648a10b22SPatrick Sanan   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetSubcommType_C",PCTelescopeSetSubcommType_Telescope);CHKERRQ(ierr);
10771e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetReductionFactor_C",PCTelescopeGetReductionFactor_Telescope);CHKERRQ(ierr);
10781e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetReductionFactor_C",PCTelescopeSetReductionFactor_Telescope);CHKERRQ(ierr);
10791e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreDM_C",PCTelescopeGetIgnoreDM_Telescope);CHKERRQ(ierr);
10801e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreDM_C",PCTelescopeSetIgnoreDM_Telescope);CHKERRQ(ierr);
10810ae7c45bSDave May   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreKSPComputeOperators_C",PCTelescopeGetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr);
10820ae7c45bSDave May   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreKSPComputeOperators_C",PCTelescopeSetIgnoreKSPComputeOperators_Telescope);CHKERRQ(ierr);
10831e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetDM_C",PCTelescopeGetDM_Telescope);CHKERRQ(ierr);
10841e07b27eSBarry Smith   PetscFunctionReturn(0);
10851e07b27eSBarry Smith }
1086