xref: /petsc/src/ksp/pc/impls/telescope/telescope.c (revision bd49479cc7f900be4c8295459593dede491c20c2)
11e07b27eSBarry Smith 
21e07b27eSBarry Smith /*
31e07b27eSBarry Smith 
41e07b27eSBarry Smith  Defines a telescoping preconditioner
51e07b27eSBarry Smith 
61e07b27eSBarry Smith  Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation
71e07b27eSBarry Smith  creates a child sub-communicator (c') containing less MPI processes than the original parent preconditioner (PC).
81e07b27eSBarry Smith 
91e07b27eSBarry Smith  During PCSetup, the B operator is scattered onto c'.
101e07b27eSBarry Smith  Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c').
111e07b27eSBarry Smith  Then KSPSolve() is executed on the c' communicator.
121e07b27eSBarry Smith 
131e07b27eSBarry Smith  The preconditioner is deemed telescopint as it only calls KSPSolve() on a single
141e07b27eSBarry Smith  sub-communicator in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators.
151e07b27eSBarry Smith  This means there will be MPI processes within c, which will be idle during the application of this preconditioner.
161e07b27eSBarry Smith 
171e07b27eSBarry Smith  Comments:
181e07b27eSBarry Smith  - The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED
191e07b27eSBarry Smith  creation routine. We run the sub KSP on only the ranks within the communicator which have a color equal to zero.
201e07b27eSBarry Smith 
211e07b27eSBarry Smith  - The telescoping preconditioner is aware of nullspaces which are attached to the only B operator.
221e07b27eSBarry Smith  In case where B has a n nullspace attached, these nullspaces vectors are extract from B and mapped into
231e07b27eSBarry Smith  a new nullspace (defined on the sub-communicator) which is attached to B' (the B operator which was scattered to c')
241e07b27eSBarry Smith 
251e07b27eSBarry Smith  - The telescoping preconditioner is aware of an attached DM. In the event that the DM is of type DMDA (2D or 3D -
261e07b27eSBarry 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
271e07b27eSBarry Smith  is attached the sub KSPSolve(). The design of telescope is such that it should be possible to extend support
281e07b27eSBarry Smith  for re-partitioning other DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs.
291e07b27eSBarry Smith 
301e07b27eSBarry Smith  - By default, B' is defined by simply fusing rows from different MPI processes
311e07b27eSBarry Smith 
321e07b27eSBarry Smith  - When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permuting of B
331e07b27eSBarry Smith  into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatGetSubMatrices(), (iii) fusing the
341e07b27eSBarry Smith  locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat()
351e07b27eSBarry Smith 
361e07b27eSBarry Smith  Limitations/improvements
371e07b27eSBarry Smith  - VecPlaceArray could be used within PCApply() to improve efficiency and reduce memory usage.
381e07b27eSBarry Smith 
391e07b27eSBarry Smith  - The symmetric permutation used when a DMDA is encountered is performed via explicitly assmbleming a permutation matrix P,
401e07b27eSBarry 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
411e07b27eSBarry Smith  VecPermute() does not supported for the use case required here. By computing P, I can permute both the operator and RHS in a
421e07b27eSBarry Smith  consistent manner.
431e07b27eSBarry Smith 
441e07b27eSBarry Smith  - Mapping of vectors is performed this way
451e07b27eSBarry 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.
461e07b27eSBarry Smith  Using the interlaced creation routine, the ranks in c with color = 0, will be rank 0 and 2.
471e07b27eSBarry Smith  We perform the scatter to the sub-comm in the following way,
481e07b27eSBarry Smith  [1] Given a vector x defined on comm c
491e07b27eSBarry Smith 
501e07b27eSBarry Smith    rank(c) : _________ 0 ______  ________ 1 _______  ________ 2 _____________ ___________ 3 __________
511e07b27eSBarry 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]
521e07b27eSBarry Smith 
531e07b27eSBarry Smith  scatter to xtmp defined also on comm c so that we have the following values
541e07b27eSBarry Smith 
551e07b27eSBarry Smith    rank(c) : ___________________ 0 ________________  _1_  ______________________ 2 _______________________  __3_
561e07b27eSBarry Smith       xtmp : [0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11] [ 6 ] [12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23] [ 18 ]
571e07b27eSBarry Smith 
581e07b27eSBarry Smith  The entries on rank 1 and 3 (ranks which do not have a color = 0 in c') are filled with the first bs values living on that rank,
591e07b27eSBarry Smith  where bs is the block size of the vector.
601e07b27eSBarry Smith 
611e07b27eSBarry 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'.
621e07b27eSBarry Smith  Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0.
631e07b27eSBarry Smith 
641e07b27eSBarry Smith   rank(c') : ___________________ 0 _______________  ______________________ 1 _____________________
651e07b27eSBarry 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]
661e07b27eSBarry Smith 
671e07b27eSBarry Smith */
681e07b27eSBarry Smith 
691e07b27eSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
701e07b27eSBarry Smith #include <petscksp.h> /*I "petscksp.h" I*/
711e07b27eSBarry Smith #include <petscdm.h> /*I "petscdm.h" I*/
721e07b27eSBarry Smith 
731e07b27eSBarry Smith #include "telescope.h"
741e07b27eSBarry Smith 
751e07b27eSBarry Smith /*
761e07b27eSBarry Smith  PCTelescopeSetUp_default()
771e07b27eSBarry Smith  PCTelescopeMatCreate_default()
781e07b27eSBarry Smith 
791e07b27eSBarry Smith  default
801e07b27eSBarry Smith 
811e07b27eSBarry Smith  // scatter in
821e07b27eSBarry Smith  x(comm) -> xtmp(comm)
831e07b27eSBarry Smith 
841e07b27eSBarry Smith  xred(subcomm) <- xtmp
851e07b27eSBarry Smith  yred(subcomm)
861e07b27eSBarry Smith 
871e07b27eSBarry Smith  yred(subcomm) --> xtmp
881e07b27eSBarry Smith 
891e07b27eSBarry Smith  // scatter out
901e07b27eSBarry Smith  xtmp(comm) -> y(comm)
911e07b27eSBarry Smith */
921e07b27eSBarry Smith 
931e07b27eSBarry Smith PetscBool isActiveRank(PetscSubcomm scomm)
941e07b27eSBarry Smith {
951e07b27eSBarry Smith   if (scomm->color == 0) { return PETSC_TRUE; }
961e07b27eSBarry Smith   else { return PETSC_FALSE; }
971e07b27eSBarry Smith }
981e07b27eSBarry Smith 
991e07b27eSBarry Smith #undef __FUNCT__
1001e07b27eSBarry Smith #define __FUNCT__ "private_PCTelescopeGetSubDM"
1011e07b27eSBarry Smith DM private_PCTelescopeGetSubDM(PC_Telescope sred)
1021e07b27eSBarry Smith {
103c6a0d831SBarry Smith   DM subdm = NULL;
1041e07b27eSBarry Smith 
1051e07b27eSBarry Smith   if (!isActiveRank(sred->psubcomm)) { subdm = NULL; }
1061e07b27eSBarry Smith   else {
1071e07b27eSBarry Smith     switch (sred->sr_type) {
1081e07b27eSBarry Smith     case TELESCOPE_DEFAULT: subdm = NULL;
1091e07b27eSBarry Smith       break;
1101e07b27eSBarry Smith     case TELESCOPE_DMDA:    subdm = ((PC_Telescope_DMDACtx*)sred->dm_ctx)->dmrepart;
1111e07b27eSBarry Smith       break;
1121e07b27eSBarry Smith     case TELESCOPE_DMPLEX:  subdm = NULL;
1131e07b27eSBarry Smith       break;
1141e07b27eSBarry Smith     }
1151e07b27eSBarry Smith   }
1161e07b27eSBarry Smith   return(subdm);
1171e07b27eSBarry Smith }
1181e07b27eSBarry Smith 
1191e07b27eSBarry Smith #undef __FUNCT__
1201e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeSetUp_default"
1211e07b27eSBarry Smith PetscErrorCode PCTelescopeSetUp_default(PC pc,PC_Telescope sred)
1221e07b27eSBarry Smith {
1231e07b27eSBarry Smith   PetscErrorCode ierr;
1241e07b27eSBarry Smith   PetscInt       m,M,bs,st,ed;
1251e07b27eSBarry Smith   Vec            x,xred,yred,xtmp;
1261e07b27eSBarry Smith   Mat            B;
1271e07b27eSBarry Smith   MPI_Comm       comm,subcomm;
1281e07b27eSBarry Smith   VecScatter     scatter;
1291e07b27eSBarry Smith   IS             isin;
1301e07b27eSBarry Smith 
1311e07b27eSBarry Smith   PetscFunctionBegin;
1321e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: setup (default)\n");CHKERRQ(ierr);
1331e07b27eSBarry Smith   comm = PetscSubcommParent(sred->psubcomm);
1341e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
1351e07b27eSBarry Smith 
1361e07b27eSBarry Smith   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
1371e07b27eSBarry Smith   ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr);
1381e07b27eSBarry Smith   ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr);
1391e07b27eSBarry Smith   ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr);
1401e07b27eSBarry Smith 
1411e07b27eSBarry Smith   xred = NULL;
1421e07b27eSBarry Smith   m = bs;
1431e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
1441e07b27eSBarry Smith     ierr = VecCreate(subcomm,&xred);CHKERRQ(ierr);
1451e07b27eSBarry Smith     ierr = VecSetSizes(xred,PETSC_DECIDE,M);CHKERRQ(ierr);
1461e07b27eSBarry Smith     ierr = VecSetBlockSize(xred,bs);CHKERRQ(ierr);
1471e07b27eSBarry Smith     ierr = VecSetFromOptions(xred);CHKERRQ(ierr);
1481e07b27eSBarry Smith     ierr = VecGetLocalSize(xred,&m);
1491e07b27eSBarry Smith   }
1501e07b27eSBarry Smith 
1511e07b27eSBarry Smith   yred = NULL;
1521e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
1531e07b27eSBarry Smith     ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr);
1541e07b27eSBarry Smith   }
1551e07b27eSBarry Smith 
1561e07b27eSBarry Smith   ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr);
1571e07b27eSBarry Smith   ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr);
1581e07b27eSBarry Smith   ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr);
1591e07b27eSBarry Smith   ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr);
1601e07b27eSBarry Smith 
1611e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
1621e07b27eSBarry Smith     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
1631e07b27eSBarry Smith     ierr = ISCreateStride(comm,(ed-st),st,1,&isin);CHKERRQ(ierr);
1641e07b27eSBarry Smith   } else {
1651e07b27eSBarry Smith     ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr);
1661e07b27eSBarry Smith     ierr = ISCreateStride(comm,bs,st,1,&isin);CHKERRQ(ierr);
1671e07b27eSBarry Smith   }
1681e07b27eSBarry Smith   ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr);
1691e07b27eSBarry Smith 
1701e07b27eSBarry Smith   ierr = VecScatterCreate(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr);
1711e07b27eSBarry Smith 
1721e07b27eSBarry Smith   sred->isin    = isin;
1731e07b27eSBarry Smith   sred->scatter = scatter;
1741e07b27eSBarry Smith   sred->xred    = xred;
1751e07b27eSBarry Smith   sred->yred    = yred;
1761e07b27eSBarry Smith   sred->xtmp    = xtmp;
1771e07b27eSBarry Smith   ierr = VecDestroy(&x);CHKERRQ(ierr);
1781e07b27eSBarry Smith   PetscFunctionReturn(0);
1791e07b27eSBarry Smith }
1801e07b27eSBarry Smith 
1811e07b27eSBarry Smith #undef __FUNCT__
1821e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeMatCreate_default"
1831e07b27eSBarry Smith PetscErrorCode PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A)
1841e07b27eSBarry Smith {
1851e07b27eSBarry Smith   PetscErrorCode ierr;
1861e07b27eSBarry Smith   MPI_Comm       comm,subcomm;
1871e07b27eSBarry Smith   Mat            Bred,B;
1881e07b27eSBarry Smith   PetscInt       nr,nc;
1891e07b27eSBarry Smith   IS             isrow,iscol;
1901e07b27eSBarry Smith   Mat            Blocal,*_Blocal;
1911e07b27eSBarry Smith 
1921e07b27eSBarry Smith   PetscFunctionBegin;
1931e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (default)\n");CHKERRQ(ierr);
1941e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1951e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
1961e07b27eSBarry Smith   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
1971e07b27eSBarry Smith   ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr);
1981e07b27eSBarry Smith   isrow = sred->isin;
1991e07b27eSBarry Smith   ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr);
2001e07b27eSBarry Smith   ierr = MatGetSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr);
2011e07b27eSBarry Smith   Blocal = *_Blocal;
2021e07b27eSBarry Smith   ierr = PetscFree(_Blocal);CHKERRQ(ierr);
2031e07b27eSBarry Smith   Bred = NULL;
2041e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
2051e07b27eSBarry Smith     PetscInt mm;
2061e07b27eSBarry Smith 
2071e07b27eSBarry Smith     if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; }
2081e07b27eSBarry Smith 
2091e07b27eSBarry Smith     ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr);
2101e07b27eSBarry Smith     ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr);
2111e07b27eSBarry Smith   }
2121e07b27eSBarry Smith   *A = Bred;
2131e07b27eSBarry Smith   ierr = ISDestroy(&iscol);CHKERRQ(ierr);
2141e07b27eSBarry Smith   ierr = MatDestroy(&Blocal);CHKERRQ(ierr);
2151e07b27eSBarry Smith   PetscFunctionReturn(0);
2161e07b27eSBarry Smith }
2171e07b27eSBarry Smith 
2181e07b27eSBarry Smith #undef __FUNCT__
2191e07b27eSBarry Smith #define __FUNCT__ "PCTelescopeMatNullSpaceCreate_default"
2201e07b27eSBarry Smith PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat)
2211e07b27eSBarry Smith {
2221e07b27eSBarry Smith   PetscErrorCode   ierr;
2231e07b27eSBarry Smith   MatNullSpace     nullspace,sub_nullspace;
2241e07b27eSBarry Smith   Mat              A,B;
2251e07b27eSBarry Smith   PetscBool        has_const;
2261e07b27eSBarry Smith   PetscInt         i,k,n;
2271e07b27eSBarry Smith   const Vec        *vecs;
2281e07b27eSBarry Smith   Vec              *sub_vecs;
2291e07b27eSBarry Smith   MPI_Comm         subcomm;
2301e07b27eSBarry Smith 
2311e07b27eSBarry Smith   PetscFunctionBegin;
2321e07b27eSBarry Smith   ierr = PCGetOperators(pc,&A,&B);CHKERRQ(ierr);
2331e07b27eSBarry Smith   ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr);
2341e07b27eSBarry Smith   if (!nullspace) return(0);
2351e07b27eSBarry Smith 
2361e07b27eSBarry Smith   ierr = PetscInfo(pc,"PCTelescope: generating nullspace (default)\n");CHKERRQ(ierr);
2371e07b27eSBarry Smith   subcomm = PetscSubcommChild(sred->psubcomm);
2381e07b27eSBarry Smith   ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr);
2391e07b27eSBarry Smith 
2401e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
2411e07b27eSBarry Smith     sub_vecs = NULL;
2421e07b27eSBarry Smith     /* create new vectors */
2431e07b27eSBarry Smith     if (n != 0) {
2441e07b27eSBarry Smith       PetscMalloc(sizeof(Vec)*n,&sub_vecs);
2451e07b27eSBarry Smith       for (k=0; k<n; k++) {
2461e07b27eSBarry Smith         ierr = VecDuplicate(sred->xred,&sub_vecs[k]);CHKERRQ(ierr);
2471e07b27eSBarry Smith       }
2481e07b27eSBarry Smith     }
2491e07b27eSBarry Smith   }
2501e07b27eSBarry Smith 
2511e07b27eSBarry Smith   /* copy entries */
2521e07b27eSBarry Smith   for (k=0; k<n; k++) {
2531e07b27eSBarry Smith     const PetscScalar *x_array;
2541e07b27eSBarry Smith     PetscScalar       *LA_sub_vec;
2551e07b27eSBarry Smith     PetscInt          st,ed,bs;
2561e07b27eSBarry Smith 
2571e07b27eSBarry Smith     /* pull in vector x->xtmp */
2581e07b27eSBarry Smith     ierr = VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2591e07b27eSBarry Smith     ierr = VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2601e07b27eSBarry Smith     /* copy vector entires into xred */
2611e07b27eSBarry Smith     ierr = VecGetBlockSize(sred->xtmp,&bs);CHKERRQ(ierr);
2621e07b27eSBarry Smith     ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
2631e07b27eSBarry Smith     if (sub_vecs[k]) {
2641e07b27eSBarry Smith       ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr);
2651e07b27eSBarry Smith       ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
2661e07b27eSBarry Smith       for (i=0; i<ed-st; i++) {
2671e07b27eSBarry Smith         LA_sub_vec[i] = x_array[i];
2681e07b27eSBarry Smith       }
2691e07b27eSBarry Smith       ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
2701e07b27eSBarry Smith     }
2711e07b27eSBarry Smith     ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
2721e07b27eSBarry Smith   }
2731e07b27eSBarry Smith 
2741e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
2751e07b27eSBarry Smith     /* create new nullspace for redundant object */
2761e07b27eSBarry Smith     ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,&sub_nullspace);CHKERRQ(ierr);
2771e07b27eSBarry Smith     /* attach redundant nullspace to Bred */
2781e07b27eSBarry Smith     ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr);
2791e07b27eSBarry Smith 
2801e07b27eSBarry Smith     for (k=0; k<n; k++) {
2811e07b27eSBarry Smith       ierr = VecDestroy(&sub_vecs[k]);CHKERRQ(ierr);
2821e07b27eSBarry Smith     }
2831e07b27eSBarry Smith     if (sub_vecs) PetscFree(sub_vecs);
2841e07b27eSBarry Smith   }
2851e07b27eSBarry Smith   PetscFunctionReturn(0);
2861e07b27eSBarry Smith }
2871e07b27eSBarry Smith 
2881e07b27eSBarry Smith #undef __FUNCT__
2891e07b27eSBarry Smith #define __FUNCT__ "PCView_Telescope"
2901e07b27eSBarry Smith static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer)
2911e07b27eSBarry Smith {
2921e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
2931e07b27eSBarry Smith   PetscErrorCode   ierr;
2941e07b27eSBarry Smith   PetscBool        iascii,isstring;
2951e07b27eSBarry Smith   PetscViewer      subviewer;
2961e07b27eSBarry Smith 
2971e07b27eSBarry Smith   PetscFunctionBegin;
2981e07b27eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
2991e07b27eSBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
3001e07b27eSBarry Smith   if (iascii) {
3011e07b27eSBarry Smith     if (!sred->psubcomm) {
3021e07b27eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: preconditioner not yet setup\n");CHKERRQ(ierr);
3031e07b27eSBarry Smith     } else {
3041e07b27eSBarry Smith       MPI_Comm    comm,subcomm;
3051e07b27eSBarry Smith       PetscMPIInt comm_size,subcomm_size;
3061e07b27eSBarry Smith       DM          dm,subdm;
3071e07b27eSBarry Smith 
3081e07b27eSBarry Smith       ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
3091e07b27eSBarry Smith       subdm = private_PCTelescopeGetSubDM(sred);
3101e07b27eSBarry Smith       comm = PetscSubcommParent(sred->psubcomm);
3111e07b27eSBarry Smith       subcomm = PetscSubcommChild(sred->psubcomm);
3121e07b27eSBarry Smith       ierr = MPI_Comm_size(comm,&comm_size);CHKERRQ(ierr);
3131e07b27eSBarry Smith       ierr = MPI_Comm_size(subcomm,&subcomm_size);CHKERRQ(ierr);
3141e07b27eSBarry Smith 
3151e07b27eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: parent comm size reduction factor = %D\n",sred->redfactor);CHKERRQ(ierr);
3161e07b27eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: comm_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size);CHKERRQ(ierr);
3171e07b27eSBarry Smith       ierr = PetscViewerGetSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
3181e07b27eSBarry Smith       if (isActiveRank(sred->psubcomm)) {
3191e07b27eSBarry Smith         ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr);
3201e07b27eSBarry Smith 
3211e07b27eSBarry Smith         if (dm && sred->ignore_dm) {
3221e07b27eSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: ignoring DM\n");CHKERRQ(ierr);
3231e07b27eSBarry Smith         }
3241e07b27eSBarry Smith         switch (sred->sr_type) {
3251e07b27eSBarry Smith         case TELESCOPE_DEFAULT:
3261e07b27eSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: using default setup\n");CHKERRQ(ierr);
3271e07b27eSBarry Smith           break;
3281e07b27eSBarry Smith         case TELESCOPE_DMDA:
3291e07b27eSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: DMDA detected\n");CHKERRQ(ierr);
3301e07b27eSBarry Smith           ierr = DMView_DMDAShort(subdm,subviewer);CHKERRQ(ierr);
3311e07b27eSBarry Smith           break;
3321e07b27eSBarry Smith         case TELESCOPE_DMPLEX:
3331e07b27eSBarry Smith           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: DMPLEX detected\n");CHKERRQ(ierr);
3341e07b27eSBarry Smith           break;
3351e07b27eSBarry Smith         }
3361e07b27eSBarry Smith 
3371e07b27eSBarry Smith         ierr = KSPView(sred->ksp,subviewer);CHKERRQ(ierr);
3381e07b27eSBarry Smith         ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr);
3391e07b27eSBarry Smith       }
3401e07b27eSBarry Smith       ierr = PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
3411e07b27eSBarry Smith     }
3421e07b27eSBarry Smith   }
3431e07b27eSBarry Smith   PetscFunctionReturn(0);
3441e07b27eSBarry Smith }
3451e07b27eSBarry Smith 
3461e07b27eSBarry Smith #undef __FUNCT__
3471e07b27eSBarry Smith #define __FUNCT__ "PCSetUp_Telescope"
3481e07b27eSBarry Smith static PetscErrorCode PCSetUp_Telescope(PC pc)
3491e07b27eSBarry Smith {
3501e07b27eSBarry Smith   PC_Telescope      sred = (PC_Telescope)pc->data;
3511e07b27eSBarry Smith   PetscErrorCode    ierr;
352*bd49479cSSatish Balay   MPI_Comm          comm,subcomm=0;
3531e07b27eSBarry Smith   PCTelescopeType   sr_type;
3541e07b27eSBarry Smith 
3551e07b27eSBarry Smith   PetscFunctionBegin;
3561e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
3571e07b27eSBarry Smith 
3581e07b27eSBarry Smith   /* subcomm definition */
3591e07b27eSBarry Smith   if (!pc->setupcalled) {
3601e07b27eSBarry Smith     if (!sred->psubcomm) {
3611e07b27eSBarry Smith       ierr = PetscSubcommCreate(comm,&sred->psubcomm);CHKERRQ(ierr);
3621e07b27eSBarry Smith       ierr = PetscSubcommSetNumber(sred->psubcomm,sred->redfactor);CHKERRQ(ierr);
3631e07b27eSBarry Smith       ierr = PetscSubcommSetType(sred->psubcomm,PETSC_SUBCOMM_INTERLACED);CHKERRQ(ierr);
3641e07b27eSBarry Smith       /* disable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */
3651e07b27eSBarry Smith       /* ierr = PetscSubcommSetFromOptions(sred->psubcomm);CHKERRQ(ierr); */
3661e07b27eSBarry Smith       ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr);
3671e07b27eSBarry Smith       /* create a new PC that processors in each subcomm have copy of */
3681e07b27eSBarry Smith     }
3691e07b27eSBarry Smith   }
3700f6f40a7SSatish Balay   subcomm = PetscSubcommChild(sred->psubcomm);
3711e07b27eSBarry Smith 
3721e07b27eSBarry Smith   /* internal KSP */
3731e07b27eSBarry Smith   if (!pc->setupcalled) {
3741e07b27eSBarry Smith     const char *prefix;
3751e07b27eSBarry Smith 
3761e07b27eSBarry Smith     if (isActiveRank(sred->psubcomm)) {
3771e07b27eSBarry Smith       ierr = KSPCreate(subcomm,&sred->ksp);CHKERRQ(ierr);
3781e07b27eSBarry Smith       ierr = KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure);CHKERRQ(ierr);
3791e07b27eSBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
3801e07b27eSBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp);CHKERRQ(ierr);
3811e07b27eSBarry Smith       ierr = KSPSetType(sred->ksp,KSPPREONLY);CHKERRQ(ierr);
3821e07b27eSBarry Smith       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
3831e07b27eSBarry Smith       ierr = KSPSetOptionsPrefix(sred->ksp,prefix);CHKERRQ(ierr);
3841e07b27eSBarry Smith       ierr = KSPAppendOptionsPrefix(sred->ksp,"telescope_");CHKERRQ(ierr);
3851e07b27eSBarry Smith     }
3861e07b27eSBarry Smith   }
3871e07b27eSBarry Smith   /* Determine type of setup/update */
3881e07b27eSBarry Smith   if (!pc->setupcalled) {
3891e07b27eSBarry Smith     PetscBool has_dm,same;
3901e07b27eSBarry Smith     DM        dm;
3911e07b27eSBarry Smith 
3921e07b27eSBarry Smith     sr_type = TELESCOPE_DEFAULT;
3931e07b27eSBarry Smith     has_dm = PETSC_FALSE;
3941e07b27eSBarry Smith     ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
3951e07b27eSBarry Smith     if (dm) { has_dm = PETSC_TRUE; }
3961e07b27eSBarry Smith     if (has_dm) {
3971e07b27eSBarry Smith       /* check for dmda */
3981e07b27eSBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&same);CHKERRQ(ierr);
3991e07b27eSBarry Smith       if (same) {
4001e07b27eSBarry Smith         ierr = PetscInfo(pc,"PCTelescope: found DMDA\n");CHKERRQ(ierr);
4011e07b27eSBarry Smith         sr_type = TELESCOPE_DMDA;
4021e07b27eSBarry Smith       }
4031e07b27eSBarry Smith       /* check for dmplex */
4041e07b27eSBarry Smith       ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same);CHKERRQ(ierr);
4051e07b27eSBarry Smith       if (same) {
4061e07b27eSBarry Smith         PetscInfo(pc,"PCTelescope: found DMPLEX\n");
4071e07b27eSBarry Smith         sr_type = TELESCOPE_DMPLEX;
4081e07b27eSBarry Smith       }
4091e07b27eSBarry Smith     }
4101e07b27eSBarry Smith 
4111e07b27eSBarry Smith     if (sred->ignore_dm) {
4121e07b27eSBarry Smith       PetscInfo(pc,"PCTelescope: ignore DM\n");
4131e07b27eSBarry Smith       sr_type = TELESCOPE_DEFAULT;
4141e07b27eSBarry Smith     }
4151e07b27eSBarry Smith     sred->sr_type = sr_type;
4161e07b27eSBarry Smith   } else {
4171e07b27eSBarry Smith     sr_type = sred->sr_type;
4181e07b27eSBarry Smith   }
4191e07b27eSBarry Smith 
4201e07b27eSBarry Smith   /* set function pointers for repartition setup, matrix creation/update, matrix nullspace and reset functionality */
4211e07b27eSBarry Smith   switch (sr_type) {
4221e07b27eSBarry Smith   case TELESCOPE_DEFAULT:
4231e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
4241e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
4251e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
4261e07b27eSBarry Smith     sred->pctelescope_reset_type              = NULL;
4271e07b27eSBarry Smith     break;
4281e07b27eSBarry Smith   case TELESCOPE_DMDA:
4291e07b27eSBarry Smith     pc->ops->apply                          = PCApply_Telescope_dmda;
4301e07b27eSBarry Smith     sred->pctelescope_setup_type              = PCTelescopeSetUp_dmda;
4311e07b27eSBarry Smith     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_dmda;
4321e07b27eSBarry Smith     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda;
4331e07b27eSBarry Smith     sred->pctelescope_reset_type              = PCReset_Telescope_dmda;
4341e07b27eSBarry Smith     break;
4351e07b27eSBarry Smith   case TELESCOPE_DMPLEX: SETERRQ(comm,PETSC_ERR_SUP,"Supprt for DMPLEX is currently not available");
4361e07b27eSBarry Smith     break;
4371e07b27eSBarry Smith   default: SETERRQ(comm,PETSC_ERR_SUP,"Only supprt for repartitioning DMDA is provided");
4381e07b27eSBarry Smith     break;
4391e07b27eSBarry Smith   }
4401e07b27eSBarry Smith 
4411e07b27eSBarry Smith   /* setup */
4421e07b27eSBarry Smith   if (sred->pctelescope_setup_type) {
4431e07b27eSBarry Smith     ierr = sred->pctelescope_setup_type(pc,sred);CHKERRQ(ierr);
4441e07b27eSBarry Smith   }
4451e07b27eSBarry Smith   /* update */
4461e07b27eSBarry Smith   if (!pc->setupcalled) {
4471e07b27eSBarry Smith     if (sred->pctelescope_matcreate_type) {
4481e07b27eSBarry Smith       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred);CHKERRQ(ierr);
4491e07b27eSBarry Smith     }
4501e07b27eSBarry Smith     if (sred->pctelescope_matnullspacecreate_type) {
4511e07b27eSBarry Smith       ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr);
4521e07b27eSBarry Smith     }
4531e07b27eSBarry Smith   } else {
4541e07b27eSBarry Smith     if (sred->pctelescope_matcreate_type) {
4551e07b27eSBarry Smith       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred);CHKERRQ(ierr);
4561e07b27eSBarry Smith     }
4571e07b27eSBarry Smith   }
4581e07b27eSBarry Smith 
4591e07b27eSBarry Smith   /* common - no construction */
4601e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
4611e07b27eSBarry Smith     ierr = KSPSetOperators(sred->ksp,sred->Bred,sred->Bred);CHKERRQ(ierr);
4621e07b27eSBarry Smith     if (pc->setfromoptionscalled && !pc->setupcalled){
4631e07b27eSBarry Smith       ierr = KSPSetFromOptions(sred->ksp);CHKERRQ(ierr);
4641e07b27eSBarry Smith     }
4651e07b27eSBarry Smith   }
4661e07b27eSBarry Smith   PetscFunctionReturn(0);
4671e07b27eSBarry Smith }
4681e07b27eSBarry Smith 
4691e07b27eSBarry Smith #undef __FUNCT__
4701e07b27eSBarry Smith #define __FUNCT__ "PCApply_Telescope"
4711e07b27eSBarry Smith static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y)
4721e07b27eSBarry Smith {
4731e07b27eSBarry Smith   PC_Telescope      sred = (PC_Telescope)pc->data;
4741e07b27eSBarry Smith   PetscErrorCode    ierr;
4751e07b27eSBarry Smith   Vec               xtmp,xred,yred;
4761e07b27eSBarry Smith   PetscInt          i,st,ed,bs;
4771e07b27eSBarry Smith   VecScatter        scatter;
4781e07b27eSBarry Smith   PetscScalar       *array;
4791e07b27eSBarry Smith   const PetscScalar *x_array;
4801e07b27eSBarry Smith 
4811e07b27eSBarry Smith   PetscFunctionBegin;
4821e07b27eSBarry Smith   xtmp    = sred->xtmp;
4831e07b27eSBarry Smith   scatter = sred->scatter;
4841e07b27eSBarry Smith   xred    = sred->xred;
4851e07b27eSBarry Smith   yred    = sred->yred;
4861e07b27eSBarry Smith 
4871e07b27eSBarry Smith   /* pull in vector x->xtmp */
4881e07b27eSBarry Smith   ierr = VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4891e07b27eSBarry Smith   ierr = VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
4901e07b27eSBarry Smith 
4911e07b27eSBarry Smith   /* copy vector entires into xred */
4921e07b27eSBarry Smith   ierr = VecGetBlockSize(xtmp,&bs);CHKERRQ(ierr);
4931e07b27eSBarry Smith   ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr);
4941e07b27eSBarry Smith   if (xred) {
4951e07b27eSBarry Smith     PetscScalar *LA_xred;
4961e07b27eSBarry Smith     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
4971e07b27eSBarry Smith     ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr);
4981e07b27eSBarry Smith     for (i=0; i<ed-st; i++) {
4991e07b27eSBarry Smith       LA_xred[i] = x_array[i];
5001e07b27eSBarry Smith     }
5011e07b27eSBarry Smith     ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr);
5021e07b27eSBarry Smith   }
5031e07b27eSBarry Smith   ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr);
5041e07b27eSBarry Smith   /* solve */
5051e07b27eSBarry Smith   if (isActiveRank(sred->psubcomm)) {
5061e07b27eSBarry Smith     ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr);
5071e07b27eSBarry Smith   }
5081e07b27eSBarry Smith   /* return vector */
5091e07b27eSBarry Smith   ierr = VecGetBlockSize(xtmp,&bs);CHKERRQ(ierr);
5101e07b27eSBarry Smith   ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr);
5111e07b27eSBarry Smith   if (yred) {
5121e07b27eSBarry Smith     const PetscScalar *LA_yred;
5131e07b27eSBarry Smith     ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr);
5141e07b27eSBarry Smith     ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr);
5151e07b27eSBarry Smith     for (i=0; i<ed-st; i++) {
5161e07b27eSBarry Smith       array[i] = LA_yred[i];
5171e07b27eSBarry Smith     }
5181e07b27eSBarry Smith     ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr);
5191e07b27eSBarry Smith   } else {
5201e07b27eSBarry Smith     for (i=0; i<bs; i++) {
5211e07b27eSBarry Smith       array[i] = 0.0;
5221e07b27eSBarry Smith     }
5231e07b27eSBarry Smith   }
5241e07b27eSBarry Smith   ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr);
5251e07b27eSBarry Smith   ierr = VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5261e07b27eSBarry Smith   ierr = VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5271e07b27eSBarry Smith   PetscFunctionReturn(0);
5281e07b27eSBarry Smith }
5291e07b27eSBarry Smith 
5301e07b27eSBarry Smith #undef __FUNCT__
5311e07b27eSBarry Smith #define __FUNCT__ "PCReset_Telescope"
5321e07b27eSBarry Smith static PetscErrorCode PCReset_Telescope(PC pc)
5331e07b27eSBarry Smith {
5341e07b27eSBarry Smith   PC_Telescope   sred = (PC_Telescope)pc->data;
5351e07b27eSBarry Smith   PetscErrorCode ierr;
5361e07b27eSBarry Smith 
5371e07b27eSBarry Smith   ierr = ISDestroy(&sred->isin);CHKERRQ(ierr);
5381e07b27eSBarry Smith   ierr = VecScatterDestroy(&sred->scatter);CHKERRQ(ierr);
5391e07b27eSBarry Smith   if (sred->xred) { ierr = VecDestroy(&sred->xred);CHKERRQ(ierr); }
5401e07b27eSBarry Smith   if (sred->yred) { ierr = VecDestroy(&sred->yred);CHKERRQ(ierr); }
5411e07b27eSBarry Smith   if (sred->xtmp) { ierr = VecDestroy(&sred->xtmp);CHKERRQ(ierr); }
5421e07b27eSBarry Smith   if (sred->Bred) { ierr = MatDestroy(&sred->Bred);CHKERRQ(ierr); }
5431e07b27eSBarry Smith   if (sred->ksp) { ierr = KSPReset(sred->ksp);CHKERRQ(ierr); }
5441e07b27eSBarry Smith   if (sred->pctelescope_reset_type) {
5451e07b27eSBarry Smith     ierr = sred->pctelescope_reset_type(pc);CHKERRQ(ierr);
5461e07b27eSBarry Smith   }
5471e07b27eSBarry Smith   PetscFunctionReturn(0);
5481e07b27eSBarry Smith }
5491e07b27eSBarry Smith 
5501e07b27eSBarry Smith #undef __FUNCT__
5511e07b27eSBarry Smith #define __FUNCT__ "PCDestroy_Telescope"
5521e07b27eSBarry Smith static PetscErrorCode PCDestroy_Telescope(PC pc)
5531e07b27eSBarry Smith {
5541e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
5551e07b27eSBarry Smith   PetscErrorCode   ierr;
5561e07b27eSBarry Smith 
5571e07b27eSBarry Smith   PetscFunctionBegin;
5581e07b27eSBarry Smith   ierr = PCReset_Telescope(pc);CHKERRQ(ierr);
5591e07b27eSBarry Smith   if (sred->ksp) { ierr = KSPDestroy(&sred->ksp);CHKERRQ(ierr); }
5601e07b27eSBarry Smith   ierr = PetscSubcommDestroy(&sred->psubcomm);CHKERRQ(ierr);
5611e07b27eSBarry Smith   if (sred->dm_ctx) PetscFree(sred->dm_ctx);
5621e07b27eSBarry Smith   PetscFree(pc->data);
5631e07b27eSBarry Smith   PetscFunctionReturn(0);
5641e07b27eSBarry Smith }
5651e07b27eSBarry Smith 
5661e07b27eSBarry Smith #undef __FUNCT__
5671e07b27eSBarry Smith #define __FUNCT__ "PCSetFromOptions_Telescope"
5681e07b27eSBarry Smith static PetscErrorCode PCSetFromOptions_Telescope(PetscOptions *PetscOptionsObject,PC pc)
5691e07b27eSBarry Smith {
5701e07b27eSBarry Smith   PC_Telescope     sred = (PC_Telescope)pc->data;
5711e07b27eSBarry Smith   PetscErrorCode   ierr;
5721e07b27eSBarry Smith   MPI_Comm         comm;
5731e07b27eSBarry Smith   PetscMPIInt      size;
5741e07b27eSBarry Smith 
5751e07b27eSBarry Smith   PetscFunctionBegin;
5761e07b27eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
5771e07b27eSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
5781e07b27eSBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Telescope options");CHKERRQ(ierr);
5791e07b27eSBarry Smith   ierr = PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,0);CHKERRQ(ierr);
5801e07b27eSBarry Smith   if (sred->redfactor > size) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size");
5811e07b27eSBarry Smith   ierr = PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,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 
5881e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp)
5891e07b27eSBarry Smith {
5901e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
591*bd49479cSSatish Balay   PetscFunctionBegin;
5921e07b27eSBarry Smith   if (ksp) *ksp = red->ksp;
593*bd49479cSSatish Balay   PetscFunctionReturn(0);
5941e07b27eSBarry Smith }
5951e07b27eSBarry Smith 
5961e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact)
5971e07b27eSBarry Smith {
5981e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
599*bd49479cSSatish Balay   PetscFunctionBegin;
6001e07b27eSBarry Smith   if (fact) *fact = red->redfactor;
601*bd49479cSSatish Balay   PetscFunctionReturn(0);
6021e07b27eSBarry Smith }
6031e07b27eSBarry Smith 
6041e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact)
6051e07b27eSBarry Smith {
6061e07b27eSBarry Smith   PC_Telescope     red = (PC_Telescope)pc->data;
6071e07b27eSBarry Smith   PetscMPIInt      size;
6081e07b27eSBarry Smith   PetscErrorCode   ierr;
6091e07b27eSBarry Smith 
610*bd49479cSSatish Balay   PetscFunctionBegin;
6111e07b27eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
6121e07b27eSBarry Smith   if (fact <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact);
6131e07b27eSBarry Smith   if (fact > size) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact);
6141e07b27eSBarry Smith   red->redfactor = fact;
615*bd49479cSSatish Balay   PetscFunctionReturn(0);
6161e07b27eSBarry Smith }
6171e07b27eSBarry Smith 
6181e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v)
6191e07b27eSBarry Smith {
6201e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
621*bd49479cSSatish Balay   PetscFunctionBegin;
6221e07b27eSBarry Smith   if (v) *v = red->ignore_dm;
623*bd49479cSSatish Balay   PetscFunctionReturn(0);
6241e07b27eSBarry Smith }
6251e07b27eSBarry Smith static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v)
6261e07b27eSBarry Smith {
6271e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
628*bd49479cSSatish Balay   PetscFunctionBegin;
6291e07b27eSBarry Smith   red->ignore_dm = v;
630*bd49479cSSatish Balay   PetscFunctionReturn(0);
6311e07b27eSBarry Smith }
6321e07b27eSBarry Smith 
6331e07b27eSBarry Smith static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm)
6341e07b27eSBarry Smith {
6351e07b27eSBarry Smith   PC_Telescope red = (PC_Telescope)pc->data;
636*bd49479cSSatish Balay   PetscFunctionBegin;
6371e07b27eSBarry Smith   *dm = private_PCTelescopeGetSubDM(red);
638*bd49479cSSatish Balay   PetscFunctionReturn(0);
6391e07b27eSBarry Smith }
6401e07b27eSBarry Smith 
6411e07b27eSBarry Smith /*@
6421e07b27eSBarry Smith  PCTelescopeGetKSP - Gets the KSP created by the telescoping PC.
6431e07b27eSBarry Smith 
6441e07b27eSBarry Smith  Not Collective
6451e07b27eSBarry Smith 
6461e07b27eSBarry Smith  Input Parameter:
6471e07b27eSBarry Smith  .  pc - the preconditioner context
6481e07b27eSBarry Smith 
6491e07b27eSBarry Smith  Output Parameter:
6501e07b27eSBarry Smith  .  subksp - the KSP defined the smaller set of processes
6511e07b27eSBarry Smith 
6521e07b27eSBarry Smith  Level: advanced
6531e07b27eSBarry Smith 
6541e07b27eSBarry Smith  .keywords: PC, telescopting solve
6551e07b27eSBarry Smith  @*/
6561e07b27eSBarry Smith PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp)
6571e07b27eSBarry Smith {
658*bd49479cSSatish Balay   PetscErrorCode ierr;
659*bd49479cSSatish Balay   PetscFunctionBegin;
660*bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp));CHKERRQ(ierr);
661*bd49479cSSatish Balay   PetscFunctionReturn(0);
6621e07b27eSBarry Smith }
6631e07b27eSBarry Smith 
6641e07b27eSBarry Smith /*@
6651e07b27eSBarry Smith  PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by.
6661e07b27eSBarry Smith 
6671e07b27eSBarry Smith  Not Collective
6681e07b27eSBarry Smith 
6691e07b27eSBarry Smith  Input Parameter:
6701e07b27eSBarry Smith  .  pc - the preconditioner context
6711e07b27eSBarry Smith 
6721e07b27eSBarry Smith  Output Parameter:
6731e07b27eSBarry Smith  .  fact - the reduction factor
6741e07b27eSBarry Smith 
6751e07b27eSBarry Smith  Level: advanced
6761e07b27eSBarry Smith 
6771e07b27eSBarry Smith  .keywords: PC, telescoping solve
6781e07b27eSBarry Smith  @*/
6791e07b27eSBarry Smith PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact)
6801e07b27eSBarry Smith {
681*bd49479cSSatish Balay   PetscErrorCode ierr;
682*bd49479cSSatish Balay   PetscFunctionBegin;
683*bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact));CHKERRQ(ierr);
684*bd49479cSSatish Balay   PetscFunctionReturn(0);
6851e07b27eSBarry Smith }
6861e07b27eSBarry Smith 
6871e07b27eSBarry Smith /*@
6881e07b27eSBarry Smith  PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by.
6891e07b27eSBarry Smith 
6901e07b27eSBarry Smith  Not Collective
6911e07b27eSBarry Smith 
6921e07b27eSBarry Smith  Input Parameter:
6931e07b27eSBarry Smith  .  pc - the preconditioner context
6941e07b27eSBarry Smith 
6951e07b27eSBarry Smith  Output Parameter:
6961e07b27eSBarry Smith  .  fact - the reduction factor
6971e07b27eSBarry Smith 
6981e07b27eSBarry Smith  Level: advanced
6991e07b27eSBarry Smith 
7001e07b27eSBarry Smith  .keywords: PC, telescoping solve
7011e07b27eSBarry Smith  @*/
7021e07b27eSBarry Smith PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact)
7031e07b27eSBarry Smith {
704*bd49479cSSatish Balay   PetscErrorCode ierr;
705*bd49479cSSatish Balay   PetscFunctionBegin;
706*bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact));CHKERRQ(ierr);
707*bd49479cSSatish Balay   PetscFunctionReturn(0);
7081e07b27eSBarry Smith }
7091e07b27eSBarry Smith 
7101e07b27eSBarry Smith /*@
7111e07b27eSBarry Smith  PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used.
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  .  v - the flag
7201e07b27eSBarry Smith 
7211e07b27eSBarry Smith  Level: advanced
7221e07b27eSBarry Smith 
7231e07b27eSBarry Smith  .keywords: PC, telescoping solve
7241e07b27eSBarry Smith  @*/
7251e07b27eSBarry Smith PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v)
7261e07b27eSBarry Smith {
727*bd49479cSSatish Balay   PetscErrorCode ierr;
728*bd49479cSSatish Balay   PetscFunctionBegin;
729*bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v));CHKERRQ(ierr);
730*bd49479cSSatish Balay   PetscFunctionReturn(0);
7311e07b27eSBarry Smith }
7321e07b27eSBarry Smith 
7331e07b27eSBarry Smith /*@
7341e07b27eSBarry Smith  PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC.
7351e07b27eSBarry Smith 
7361e07b27eSBarry Smith  Not Collective
7371e07b27eSBarry Smith 
7381e07b27eSBarry Smith  Input Parameter:
7391e07b27eSBarry Smith  .  pc - the preconditioner context
7401e07b27eSBarry Smith 
7411e07b27eSBarry Smith  Output Parameter:
7421e07b27eSBarry Smith  .  v - Use PETSC_TRUE to ignore any DM
7431e07b27eSBarry Smith 
7441e07b27eSBarry Smith  Level: advanced
7451e07b27eSBarry Smith 
7461e07b27eSBarry Smith  .keywords: PC, telescoping solve
7471e07b27eSBarry Smith  @*/
748bfd6bcc6SSatish Balay PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v)
7491e07b27eSBarry Smith {
750*bd49479cSSatish Balay   PetscErrorCode ierr;
751*bd49479cSSatish Balay   PetscFunctionBegin;
752*bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v));CHKERRQ(ierr);
753*bd49479cSSatish Balay   PetscFunctionReturn(0);
7541e07b27eSBarry Smith }
7551e07b27eSBarry Smith 
7561e07b27eSBarry Smith /*@
7571e07b27eSBarry Smith  PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP.
7581e07b27eSBarry Smith 
7591e07b27eSBarry Smith  Not Collective
7601e07b27eSBarry Smith 
7611e07b27eSBarry Smith  Input Parameter:
7621e07b27eSBarry Smith  .  pc - the preconditioner context
7631e07b27eSBarry Smith 
7641e07b27eSBarry Smith  Output Parameter:
7651e07b27eSBarry Smith  .  subdm - The re-partitioned DM
7661e07b27eSBarry Smith 
7671e07b27eSBarry Smith  Level: advanced
7681e07b27eSBarry Smith 
7691e07b27eSBarry Smith  .keywords: PC, telescoping solve
7701e07b27eSBarry Smith  @*/
7711e07b27eSBarry Smith PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm)
7721e07b27eSBarry Smith {
773*bd49479cSSatish Balay   PetscErrorCode ierr;
774*bd49479cSSatish Balay   PetscFunctionBegin;
775*bd49479cSSatish Balay   ierr = PetscTryMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm));CHKERRQ(ierr);
776*bd49479cSSatish Balay   PetscFunctionReturn(0);
7771e07b27eSBarry Smith }
7781e07b27eSBarry Smith 
7791e07b27eSBarry Smith /* -------------------------------------------------------------------------------------*/
7801e07b27eSBarry Smith /*MC
7811e07b27eSBarry Smith    PCTELESCOPE - Runs a KSP solver on a sub-group of processors. MPI processes not in the sub-communicator are idle during the solve.
7821e07b27eSBarry Smith 
7831e07b27eSBarry Smith    Options Database:
7841e07b27eSBarry Smith +  -pc_telescope_reduction_factor <n> - factor to use communicator size by, for example if you are using 64 MPI processes and
7851e07b27eSBarry Smith    use an n of 4, the new sub-communicator will be 4 defined with 64/4 processes
7861e07b27eSBarry Smith -  -pc_telescope_ignore_dm <false> - flag to indicate whether an attached DM should be ignored
7871e07b27eSBarry Smith 
7881e07b27eSBarry Smith    Level: advanced
7891e07b27eSBarry Smith 
7901e07b27eSBarry Smith    Notes:
7911e07b27eSBarry Smith    The default KSP is PREONLY. If a DM is attached to the PC, it is re-partitioned on the sub-communicator.
7921e07b27eSBarry 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.
7931e07b27eSBarry Smith    Currently only support for re-partitioning a DMDA is provided.
7941e07b27eSBarry Smith    Any nullspace attached to the original Bmat operator are extracted, re-partitioned and set on the repartitioned Bmat operator.
7951e07b27eSBarry Smith    KSPSetComputeOperators() is not propagated to the sub KSP.
7961e07b27eSBarry Smith    Currently there is no support for the flag -pc_use_amat
7971e07b27eSBarry Smith 
7981e07b27eSBarry Smith    Optimizations:
7991e07b27eSBarry Smith    (i) VecPlaceArray() could be used for scatters between the vectors defined on different sized communicators, thereby slightly reducing memory footprint;
8001e07b27eSBarry Smith    (ii) Memory re-use and faster set-up would follow if the result of P^T.A.P was not re-allocated each time PCSetUp_Telescope() was called (DMDA).
8011e07b27eSBarry Smith 
8021e07b27eSBarry Smith  Contributed by Dave May
8031e07b27eSBarry Smith 
8041e07b27eSBarry Smith .seealso:  PCTelescopeGetKSP(), PCTelescopeGetDM(),
8051e07b27eSBarry Smith  PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(),
8061e07b27eSBarry Smith  PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT
8071e07b27eSBarry Smith M*/
8081e07b27eSBarry Smith #undef __FUNCT__
8091e07b27eSBarry Smith #define __FUNCT__ "PCCreate_Telescope"
8101e07b27eSBarry Smith PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc)
8111e07b27eSBarry Smith {
8121e07b27eSBarry Smith   PetscErrorCode       ierr;
8131e07b27eSBarry Smith   struct _PC_Telescope *sred;
8141e07b27eSBarry Smith 
8151e07b27eSBarry Smith   PetscFunctionBegin;
8161e07b27eSBarry Smith   ierr = PetscNewLog(pc,&sred);CHKERRQ(ierr);
8171e07b27eSBarry Smith   sred->redfactor      = 1;
8181e07b27eSBarry Smith   sred->ignore_dm      = PETSC_FALSE;
8191e07b27eSBarry Smith   pc->data             = (void*)sred;
8201e07b27eSBarry Smith 
8211e07b27eSBarry Smith   pc->ops->apply          = PCApply_Telescope;
8221e07b27eSBarry Smith   pc->ops->applytranspose = NULL;
8231e07b27eSBarry Smith   pc->ops->setup          = PCSetUp_Telescope;
8241e07b27eSBarry Smith   pc->ops->destroy        = PCDestroy_Telescope;
8251e07b27eSBarry Smith   pc->ops->reset          = PCReset_Telescope;
8261e07b27eSBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_Telescope;
8271e07b27eSBarry Smith   pc->ops->view           = PCView_Telescope;
8281e07b27eSBarry Smith 
8291e07b27eSBarry Smith   sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
8301e07b27eSBarry Smith   sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
8311e07b27eSBarry Smith   sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
8321e07b27eSBarry Smith   sred->pctelescope_reset_type              = NULL;
8331e07b27eSBarry Smith 
8341e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetKSP_C",PCTelescopeGetKSP_Telescope);CHKERRQ(ierr);
8351e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetReductionFactor_C",PCTelescopeGetReductionFactor_Telescope);CHKERRQ(ierr);
8361e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetReductionFactor_C",PCTelescopeSetReductionFactor_Telescope);CHKERRQ(ierr);
8371e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreDM_C",PCTelescopeGetIgnoreDM_Telescope);CHKERRQ(ierr);
8381e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreDM_C",PCTelescopeSetIgnoreDM_Telescope);CHKERRQ(ierr);
8391e07b27eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetDM_C",PCTelescopeGetDM_Telescope);CHKERRQ(ierr);
8401e07b27eSBarry Smith   PetscFunctionReturn(0);
8411e07b27eSBarry Smith }
842