xref: /petsc/src/ksp/pc/impls/telescope/telescope.c (revision bfd6bcc6ac8e2f0022fa8b00b8cd29f3a32cf260)
1 
2 /*
3 
4  Defines a telescoping preconditioner
5 
6  Assuming that the parent preconditioner (PC) is defined on a communicator c, this implementation
7  creates a child sub-communicator (c') containing less MPI processes than the original parent preconditioner (PC).
8 
9  During PCSetup, the B operator is scattered onto c'.
10  Within PCApply, the RHS vector (x) is scattered into a redundant vector, xred (defined on c').
11  Then KSPSolve() is executed on the c' communicator.
12 
13  The preconditioner is deemed telescopint as it only calls KSPSolve() on a single
14  sub-communicator in contrast with PCREDUNDANT which calls KSPSolve() on N sub-communicators.
15  This means there will be MPI processes within c, which will be idle during the application of this preconditioner.
16 
17  Comments:
18  - The communicator used within the telescoping preconditioner is defined by a PetscSubcomm using the INTERLACED
19  creation routine. We run the sub KSP on only the ranks within the communicator which have a color equal to zero.
20 
21  - The telescoping preconditioner is aware of nullspaces which are attached to the only B operator.
22  In case where B has a n nullspace attached, these nullspaces vectors are extract from B and mapped into
23  a new nullspace (defined on the sub-communicator) which is attached to B' (the B operator which was scattered to c')
24 
25  - The telescoping preconditioner is aware of an attached DM. In the event that the DM is of type DMDA (2D or 3D -
26  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
27  is attached the sub KSPSolve(). The design of telescope is such that it should be possible to extend support
28  for re-partitioning other DM's (e.g. DMPLEX). The user can supply a flag to ignore attached DMs.
29 
30  - By default, B' is defined by simply fusing rows from different MPI processes
31 
32  - When a DMDA is attached to the parent preconditioner, B' is defined by: (i) performing a symmetric permuting of B
33  into the ordering defined by the DMDA on c', (ii) extracting the local chunks via MatGetSubMatrices(), (iii) fusing the
34  locally (sequential) matrices defined on the ranks common to c and c' into B' using MatCreateMPIMatConcatenateSeqMat()
35 
36  Limitations/improvements
37  - VecPlaceArray could be used within PCApply() to improve efficiency and reduce memory usage.
38 
39  - The symmetric permutation used when a DMDA is encountered is performed via explicitly assmbleming a permutation matrix P,
40  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
41  VecPermute() does not supported for the use case required here. By computing P, I can permute both the operator and RHS in a
42  consistent manner.
43 
44  - Mapping of vectors is performed this way
45  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.
46  Using the interlaced creation routine, the ranks in c with color = 0, will be rank 0 and 2.
47  We perform the scatter to the sub-comm in the following way,
48  [1] Given a vector x defined on comm c
49 
50    rank(c) : _________ 0 ______  ________ 1 _______  ________ 2 _____________ ___________ 3 __________
51          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]
52 
53  scatter to xtmp defined also on comm c so that we have the following values
54 
55    rank(c) : ___________________ 0 ________________  _1_  ______________________ 2 _______________________  __3_
56       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 ]
57 
58  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,
59  where bs is the block size of the vector.
60 
61  [2] Copy the value from rank 0, 2 (indices with respect to comm c) into the vector xred which is defined on communicator c'.
62  Ranks 0 and 2 are the only ranks in the subcomm which have a color = 0.
63 
64   rank(c') : ___________________ 0 _______________  ______________________ 1 _____________________
65       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]
66 
67 */
68 
69 #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
70 #include <petscksp.h> /*I "petscksp.h" I*/
71 #include <petscdm.h> /*I "petscdm.h" I*/
72 
73 #include "telescope.h"
74 
75 /*
76  PCTelescopeSetUp_default()
77  PCTelescopeMatCreate_default()
78 
79  default
80 
81  // scatter in
82  x(comm) -> xtmp(comm)
83 
84  xred(subcomm) <- xtmp
85  yred(subcomm)
86 
87  yred(subcomm) --> xtmp
88 
89  // scatter out
90  xtmp(comm) -> y(comm)
91 */
92 
93 PetscBool isActiveRank(PetscSubcomm scomm)
94 {
95   if (scomm->color == 0) { return PETSC_TRUE; }
96   else { return PETSC_FALSE; }
97 }
98 
99 #undef __FUNCT__
100 #define __FUNCT__ "private_PCTelescopeGetSubDM"
101 DM private_PCTelescopeGetSubDM(PC_Telescope sred)
102 {
103   DM subdm;
104 
105   if (!isActiveRank(sred->psubcomm)) { subdm = NULL; }
106   else {
107     switch (sred->sr_type) {
108     case TELESCOPE_DEFAULT: subdm = NULL;
109       break;
110     case TELESCOPE_DMDA:    subdm = ((PC_Telescope_DMDACtx*)sred->dm_ctx)->dmrepart;
111       break;
112     case TELESCOPE_DMPLEX:  subdm = NULL;
113       break;
114     }
115   }
116   return(subdm);
117 }
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "PCTelescopeSetUp_default"
121 PetscErrorCode PCTelescopeSetUp_default(PC pc,PC_Telescope sred)
122 {
123   PetscErrorCode ierr;
124   PetscInt       m,M,bs,st,ed;
125   Vec            x,xred,yred,xtmp;
126   Mat            B;
127   MPI_Comm       comm,subcomm;
128   VecScatter     scatter;
129   IS             isin;
130 
131   PetscFunctionBegin;
132   ierr = PetscInfo(pc,"PCTelescope: setup (default)\n");CHKERRQ(ierr);
133   comm = PetscSubcommParent(sred->psubcomm);
134   subcomm = PetscSubcommChild(sred->psubcomm);
135 
136   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
137   ierr = MatGetSize(B,&M,NULL);CHKERRQ(ierr);
138   ierr = MatGetBlockSize(B,&bs);CHKERRQ(ierr);
139   ierr = MatCreateVecs(B,&x,NULL);CHKERRQ(ierr);
140 
141   xred = NULL;
142   m = bs;
143   if (isActiveRank(sred->psubcomm)) {
144     ierr = VecCreate(subcomm,&xred);CHKERRQ(ierr);
145     ierr = VecSetSizes(xred,PETSC_DECIDE,M);CHKERRQ(ierr);
146     ierr = VecSetBlockSize(xred,bs);CHKERRQ(ierr);
147     ierr = VecSetFromOptions(xred);CHKERRQ(ierr);
148     ierr = VecGetLocalSize(xred,&m);
149   }
150 
151   yred = NULL;
152   if (isActiveRank(sred->psubcomm)) {
153     ierr = VecDuplicate(xred,&yred);CHKERRQ(ierr);
154   }
155 
156   ierr = VecCreate(comm,&xtmp);CHKERRQ(ierr);
157   ierr = VecSetSizes(xtmp,m,PETSC_DECIDE);CHKERRQ(ierr);
158   ierr = VecSetBlockSize(xtmp,bs);CHKERRQ(ierr);
159   ierr = VecSetType(xtmp,((PetscObject)x)->type_name);CHKERRQ(ierr);
160 
161   if (isActiveRank(sred->psubcomm)) {
162     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
163     ierr = ISCreateStride(comm,(ed-st),st,1,&isin);CHKERRQ(ierr);
164   } else {
165     ierr = VecGetOwnershipRange(x,&st,&ed);CHKERRQ(ierr);
166     ierr = ISCreateStride(comm,bs,st,1,&isin);CHKERRQ(ierr);
167   }
168   ierr = ISSetBlockSize(isin,bs);CHKERRQ(ierr);
169 
170   ierr = VecScatterCreate(x,isin,xtmp,NULL,&scatter);CHKERRQ(ierr);
171 
172   sred->isin    = isin;
173   sred->scatter = scatter;
174   sred->xred    = xred;
175   sred->yred    = yred;
176   sred->xtmp    = xtmp;
177   ierr = VecDestroy(&x);CHKERRQ(ierr);
178   PetscFunctionReturn(0);
179 }
180 
181 #undef __FUNCT__
182 #define __FUNCT__ "PCTelescopeMatCreate_default"
183 PetscErrorCode PCTelescopeMatCreate_default(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A)
184 {
185   PetscErrorCode ierr;
186   MPI_Comm       comm,subcomm;
187   Mat            Bred,B;
188   PetscInt       nr,nc;
189   IS             isrow,iscol;
190   Mat            Blocal,*_Blocal;
191 
192   PetscFunctionBegin;
193   ierr = PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (default)\n");CHKERRQ(ierr);
194   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
195   subcomm = PetscSubcommChild(sred->psubcomm);
196   ierr = PCGetOperators(pc,NULL,&B);CHKERRQ(ierr);
197   ierr = MatGetSize(B,&nr,&nc);CHKERRQ(ierr);
198   isrow = sred->isin;
199   ierr = ISCreateStride(comm,nc,0,1,&iscol);CHKERRQ(ierr);
200   ierr = MatGetSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);CHKERRQ(ierr);
201   Blocal = *_Blocal;
202   ierr = PetscFree(_Blocal);CHKERRQ(ierr);
203   Bred = NULL;
204   if (isActiveRank(sred->psubcomm)) {
205     PetscInt mm;
206 
207     if (reuse != MAT_INITIAL_MATRIX) { Bred = *A; }
208 
209     ierr = MatGetSize(Blocal,&mm,NULL);CHKERRQ(ierr);
210     /* ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,PETSC_DECIDE,reuse,&Bred);CHKERRQ(ierr); */
211     ierr = MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);CHKERRQ(ierr);
212   }
213   *A = Bred;
214   ierr = ISDestroy(&iscol);CHKERRQ(ierr);
215   ierr = MatDestroy(&Blocal);CHKERRQ(ierr);
216   PetscFunctionReturn(0);
217 }
218 
219 #undef __FUNCT__
220 #define __FUNCT__ "PCTelescopeMatNullSpaceCreate_default"
221 PetscErrorCode PCTelescopeMatNullSpaceCreate_default(PC pc,PC_Telescope sred,Mat sub_mat)
222 {
223   PetscErrorCode   ierr;
224   MatNullSpace     nullspace,sub_nullspace;
225   Mat              A,B;
226   PetscBool        has_const;
227   PetscInt         i,k,n;
228   const Vec        *vecs;
229   Vec              *sub_vecs;
230   MPI_Comm         subcomm;
231 
232   PetscFunctionBegin;
233   ierr = PCGetOperators(pc,&A,&B);CHKERRQ(ierr);
234   ierr = MatGetNullSpace(B,&nullspace);CHKERRQ(ierr);
235   if (!nullspace) return(0);
236 
237   ierr = PetscInfo(pc,"PCTelescope: generating nullspace (default)\n");CHKERRQ(ierr);
238   subcomm = PetscSubcommChild(sred->psubcomm);
239   ierr = MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);CHKERRQ(ierr);
240 
241   if (isActiveRank(sred->psubcomm)) {
242     sub_vecs = NULL;
243     /* create new vectors */
244     if (n != 0) {
245       PetscMalloc(sizeof(Vec)*n,&sub_vecs);
246       for (k=0; k<n; k++) {
247         ierr = VecDuplicate(sred->xred,&sub_vecs[k]);CHKERRQ(ierr);
248       }
249     }
250   }
251 
252   /* copy entries */
253   for (k=0; k<n; k++) {
254     const PetscScalar *x_array;
255     PetscScalar       *LA_sub_vec;
256     PetscInt          st,ed,bs;
257 
258     /* pull in vector x->xtmp */
259     ierr = VecScatterBegin(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
260     ierr = VecScatterEnd(sred->scatter,vecs[k],sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
261     /* copy vector entires into xred */
262     ierr = VecGetBlockSize(sred->xtmp,&bs);CHKERRQ(ierr);
263     ierr = VecGetArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
264     if (sub_vecs[k]) {
265       ierr = VecGetOwnershipRange(sub_vecs[k],&st,&ed);CHKERRQ(ierr);
266       ierr = VecGetArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
267       for (i=0; i<ed-st; i++) {
268         LA_sub_vec[i] = x_array[i];
269       }
270       ierr = VecRestoreArray(sub_vecs[k],&LA_sub_vec);CHKERRQ(ierr);
271     }
272     ierr = VecRestoreArrayRead(sred->xtmp,&x_array);CHKERRQ(ierr);
273   }
274 
275   if (isActiveRank(sred->psubcomm)) {
276     /* create new nullspace for redundant object */
277     ierr = MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,&sub_nullspace);CHKERRQ(ierr);
278     /* attach redundant nullspace to Bred */
279     ierr = MatSetNullSpace(sub_mat,sub_nullspace);CHKERRQ(ierr);
280 
281     for (k=0; k<n; k++) {
282       ierr = VecDestroy(&sub_vecs[k]);CHKERRQ(ierr);
283     }
284     if (sub_vecs) PetscFree(sub_vecs);
285   }
286   PetscFunctionReturn(0);
287 }
288 
289 #undef __FUNCT__
290 #define __FUNCT__ "PCView_Telescope"
291 static PetscErrorCode PCView_Telescope(PC pc,PetscViewer viewer)
292 {
293   PC_Telescope     sred = (PC_Telescope)pc->data;
294   PetscErrorCode   ierr;
295   PetscBool        iascii,isstring;
296   PetscViewer      subviewer;
297 
298   PetscFunctionBegin;
299   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
300   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
301   if (iascii) {
302     if (!sred->psubcomm) {
303       ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: preconditioner not yet setup\n");CHKERRQ(ierr);
304     } else {
305       MPI_Comm    comm,subcomm;
306       PetscMPIInt comm_size,subcomm_size;
307       DM          dm,subdm;
308 
309       ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
310       subdm = private_PCTelescopeGetSubDM(sred);
311       comm = PetscSubcommParent(sred->psubcomm);
312       subcomm = PetscSubcommChild(sred->psubcomm);
313       ierr = MPI_Comm_size(comm,&comm_size);CHKERRQ(ierr);
314       ierr = MPI_Comm_size(subcomm,&subcomm_size);CHKERRQ(ierr);
315 
316       ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: parent comm size reduction factor = %D\n",sred->redfactor);CHKERRQ(ierr);
317       ierr = PetscViewerASCIIPrintf(viewer,"  Telescope: comm_size = %d , subcomm_size = %d\n",(int)comm_size,(int)subcomm_size);CHKERRQ(ierr);
318       ierr = PetscViewerGetSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
319       if (isActiveRank(sred->psubcomm)) {
320         ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr);
321 
322         if (dm && sred->ignore_dm) {
323           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: ignoring DM\n");CHKERRQ(ierr);
324         }
325         switch (sred->sr_type) {
326         case TELESCOPE_DEFAULT:
327           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: using default setup\n");CHKERRQ(ierr);
328           break;
329         case TELESCOPE_DMDA:
330           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: DMDA detected\n");CHKERRQ(ierr);
331           ierr = DMView_DMDAShort(subdm,subviewer);CHKERRQ(ierr);
332           break;
333         case TELESCOPE_DMPLEX:
334           ierr = PetscViewerASCIIPrintf(subviewer,"  Telescope: DMPLEX detected\n");CHKERRQ(ierr);
335           break;
336         }
337 
338         ierr = KSPView(sred->ksp,subviewer);CHKERRQ(ierr);
339         ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr);
340       }
341       ierr = PetscViewerRestoreSubViewer(viewer,subcomm,&subviewer);CHKERRQ(ierr);
342     }
343   }
344   PetscFunctionReturn(0);
345 }
346 
347 #undef __FUNCT__
348 #define __FUNCT__ "PCSetUp_Telescope"
349 static PetscErrorCode PCSetUp_Telescope(PC pc)
350 {
351   PC_Telescope      sred = (PC_Telescope)pc->data;
352   PetscErrorCode    ierr;
353   MPI_Comm          comm,subcomm;
354   PCTelescopeType   sr_type;
355 
356   PetscFunctionBegin;
357   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
358 
359   /* subcomm definition */
360   if (!pc->setupcalled) {
361     if (!sred->psubcomm) {
362       ierr = PetscSubcommCreate(comm,&sred->psubcomm);CHKERRQ(ierr);
363       ierr = PetscSubcommSetNumber(sred->psubcomm,sred->redfactor);CHKERRQ(ierr);
364       ierr = PetscSubcommSetType(sred->psubcomm,PETSC_SUBCOMM_INTERLACED);CHKERRQ(ierr);
365       /* disable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */
366       /* ierr = PetscSubcommSetFromOptions(sred->psubcomm);CHKERRQ(ierr); */
367       ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr);
368       /* create a new PC that processors in each subcomm have copy of */
369       subcomm = PetscSubcommChild(sred->psubcomm);
370     }
371   } else {
372     subcomm = PetscSubcommChild(sred->psubcomm);
373   }
374 
375   /* internal KSP */
376   if (!pc->setupcalled) {
377     const char *prefix;
378 
379     if (isActiveRank(sred->psubcomm)) {
380       ierr = KSPCreate(subcomm,&sred->ksp);CHKERRQ(ierr);
381       ierr = KSPSetErrorIfNotConverged(sred->ksp,pc->erroriffailure);CHKERRQ(ierr);
382       ierr = PetscObjectIncrementTabLevel((PetscObject)sred->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
383       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)sred->ksp);CHKERRQ(ierr);
384       ierr = KSPSetType(sred->ksp,KSPPREONLY);CHKERRQ(ierr);
385       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
386       ierr = KSPSetOptionsPrefix(sred->ksp,prefix);CHKERRQ(ierr);
387       ierr = KSPAppendOptionsPrefix(sred->ksp,"telescope_");CHKERRQ(ierr);
388     }
389   }
390   /* Determine type of setup/update */
391   if (!pc->setupcalled) {
392     PetscBool has_dm,same;
393     DM        dm;
394 
395     sr_type = TELESCOPE_DEFAULT;
396     has_dm = PETSC_FALSE;
397     ierr = PCGetDM(pc,&dm);CHKERRQ(ierr);
398     if (dm) { has_dm = PETSC_TRUE; }
399     if (has_dm) {
400       /* check for dmda */
401       ierr = PetscObjectTypeCompare((PetscObject)dm,DMDA,&same);CHKERRQ(ierr);
402       if (same) {
403         ierr = PetscInfo(pc,"PCTelescope: found DMDA\n");CHKERRQ(ierr);
404         sr_type = TELESCOPE_DMDA;
405       }
406       /* check for dmplex */
407       ierr = PetscObjectTypeCompare((PetscObject)dm,DMPLEX,&same);CHKERRQ(ierr);
408       if (same) {
409         PetscInfo(pc,"PCTelescope: found DMPLEX\n");
410         sr_type = TELESCOPE_DMPLEX;
411       }
412     }
413 
414     if (sred->ignore_dm) {
415       PetscInfo(pc,"PCTelescope: ignore DM\n");
416       sr_type = TELESCOPE_DEFAULT;
417     }
418     sred->sr_type = sr_type;
419   } else {
420     sr_type = sred->sr_type;
421   }
422 
423   /* set function pointers for repartition setup, matrix creation/update, matrix nullspace and reset functionality */
424   switch (sr_type) {
425   case TELESCOPE_DEFAULT:
426     sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
427     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
428     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
429     sred->pctelescope_reset_type              = NULL;
430     break;
431   case TELESCOPE_DMDA:
432     pc->ops->apply                          = PCApply_Telescope_dmda;
433     sred->pctelescope_setup_type              = PCTelescopeSetUp_dmda;
434     sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_dmda;
435     sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_dmda;
436     sred->pctelescope_reset_type              = PCReset_Telescope_dmda;
437     break;
438   case TELESCOPE_DMPLEX: SETERRQ(comm,PETSC_ERR_SUP,"Supprt for DMPLEX is currently not available");
439     break;
440   default: SETERRQ(comm,PETSC_ERR_SUP,"Only supprt for repartitioning DMDA is provided");
441     break;
442   }
443 
444   /* setup */
445   if (sred->pctelescope_setup_type) {
446     ierr = sred->pctelescope_setup_type(pc,sred);CHKERRQ(ierr);
447   }
448   /* update */
449   if (!pc->setupcalled) {
450     if (sred->pctelescope_matcreate_type) {
451       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_INITIAL_MATRIX,&sred->Bred);CHKERRQ(ierr);
452     }
453     if (sred->pctelescope_matnullspacecreate_type) {
454       ierr = sred->pctelescope_matnullspacecreate_type(pc,sred,sred->Bred);CHKERRQ(ierr);
455     }
456   } else {
457     if (sred->pctelescope_matcreate_type) {
458       ierr = sred->pctelescope_matcreate_type(pc,sred,MAT_REUSE_MATRIX,&sred->Bred);CHKERRQ(ierr);
459     }
460   }
461 
462   /* common - no construction */
463   if (isActiveRank(sred->psubcomm)) {
464     ierr = KSPSetOperators(sred->ksp,sred->Bred,sred->Bred);CHKERRQ(ierr);
465     if (pc->setfromoptionscalled && !pc->setupcalled){
466       ierr = KSPSetFromOptions(sred->ksp);CHKERRQ(ierr);
467     }
468   }
469   PetscFunctionReturn(0);
470 }
471 
472 #undef __FUNCT__
473 #define __FUNCT__ "PCApply_Telescope"
474 static PetscErrorCode PCApply_Telescope(PC pc,Vec x,Vec y)
475 {
476   PC_Telescope      sred = (PC_Telescope)pc->data;
477   PetscErrorCode    ierr;
478   Vec               xtmp,xred,yred;
479   PetscInt          i,st,ed,bs;
480   VecScatter        scatter;
481   PetscScalar       *array;
482   const PetscScalar *x_array;
483 
484   PetscFunctionBegin;
485   xtmp    = sred->xtmp;
486   scatter = sred->scatter;
487   xred    = sred->xred;
488   yred    = sred->yred;
489 
490   /* pull in vector x->xtmp */
491   ierr = VecScatterBegin(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
492   ierr = VecScatterEnd(scatter,x,xtmp,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
493 
494   /* copy vector entires into xred */
495   ierr = VecGetBlockSize(xtmp,&bs);CHKERRQ(ierr);
496   ierr = VecGetArrayRead(xtmp,&x_array);CHKERRQ(ierr);
497   if (xred) {
498     PetscScalar *LA_xred;
499     ierr = VecGetOwnershipRange(xred,&st,&ed);CHKERRQ(ierr);
500     ierr = VecGetArray(xred,&LA_xred);CHKERRQ(ierr);
501     for (i=0; i<ed-st; i++) {
502       LA_xred[i] = x_array[i];
503     }
504     ierr = VecRestoreArray(xred,&LA_xred);CHKERRQ(ierr);
505   }
506   ierr = VecRestoreArrayRead(xtmp,&x_array);CHKERRQ(ierr);
507   /* solve */
508   if (isActiveRank(sred->psubcomm)) {
509     ierr = KSPSolve(sred->ksp,xred,yred);CHKERRQ(ierr);
510   }
511   /* return vector */
512   ierr = VecGetBlockSize(xtmp,&bs);CHKERRQ(ierr);
513   ierr = VecGetArray(xtmp,&array);CHKERRQ(ierr);
514   if (yred) {
515     const PetscScalar *LA_yred;
516     ierr = VecGetOwnershipRange(yred,&st,&ed);CHKERRQ(ierr);
517     ierr = VecGetArrayRead(yred,&LA_yred);CHKERRQ(ierr);
518     for (i=0; i<ed-st; i++) {
519       array[i] = LA_yred[i];
520     }
521     ierr = VecRestoreArrayRead(yred,&LA_yred);CHKERRQ(ierr);
522   } else {
523     for (i=0; i<bs; i++) {
524       array[i] = 0.0;
525     }
526   }
527   ierr = VecRestoreArray(xtmp,&array);CHKERRQ(ierr);
528   ierr = VecScatterBegin(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
529   ierr = VecScatterEnd(scatter,xtmp,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
530   PetscFunctionReturn(0);
531 }
532 
533 #undef __FUNCT__
534 #define __FUNCT__ "PCReset_Telescope"
535 static PetscErrorCode PCReset_Telescope(PC pc)
536 {
537   PC_Telescope   sred = (PC_Telescope)pc->data;
538   PetscErrorCode ierr;
539 
540   ierr = ISDestroy(&sred->isin);CHKERRQ(ierr);
541   ierr = VecScatterDestroy(&sred->scatter);CHKERRQ(ierr);
542   if (sred->xred) { ierr = VecDestroy(&sred->xred);CHKERRQ(ierr); }
543   if (sred->yred) { ierr = VecDestroy(&sred->yred);CHKERRQ(ierr); }
544   if (sred->xtmp) { ierr = VecDestroy(&sred->xtmp);CHKERRQ(ierr); }
545   if (sred->Bred) { ierr = MatDestroy(&sred->Bred);CHKERRQ(ierr); }
546   if (sred->ksp) { ierr = KSPReset(sred->ksp);CHKERRQ(ierr); }
547   if (sred->pctelescope_reset_type) {
548     ierr = sred->pctelescope_reset_type(pc);CHKERRQ(ierr);
549   }
550   PetscFunctionReturn(0);
551 }
552 
553 #undef __FUNCT__
554 #define __FUNCT__ "PCDestroy_Telescope"
555 static PetscErrorCode PCDestroy_Telescope(PC pc)
556 {
557   PC_Telescope     sred = (PC_Telescope)pc->data;
558   PetscErrorCode   ierr;
559 
560   PetscFunctionBegin;
561   ierr = PCReset_Telescope(pc);CHKERRQ(ierr);
562   if (sred->ksp) { ierr = KSPDestroy(&sred->ksp);CHKERRQ(ierr); }
563   ierr = PetscSubcommDestroy(&sred->psubcomm);CHKERRQ(ierr);
564   if (sred->dm_ctx) PetscFree(sred->dm_ctx);
565   PetscFree(pc->data);
566   PetscFunctionReturn(0);
567 }
568 
569 #undef __FUNCT__
570 #define __FUNCT__ "PCSetFromOptions_Telescope"
571 static PetscErrorCode PCSetFromOptions_Telescope(PetscOptions *PetscOptionsObject,PC pc)
572 {
573   PC_Telescope     sred = (PC_Telescope)pc->data;
574   PetscErrorCode   ierr;
575   MPI_Comm         comm;
576   PetscMPIInt      size;
577 
578   PetscFunctionBegin;
579   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
580   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
581   ierr = PetscOptionsHead(PetscOptionsObject,"Telescope options");CHKERRQ(ierr);
582   ierr = PetscOptionsInt("-pc_telescope_reduction_factor","Factor to reduce comm size by","PCTelescopeSetReductionFactor",sred->redfactor,&sred->redfactor,0);CHKERRQ(ierr);
583   if (sred->redfactor > size) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"-pc_telescope_reduction_factor <= comm size");
584   ierr = PetscOptionsBool("-pc_telescope_ignore_dm","Ignore any DM attached to the PC","PCTelescopeSetIgnoreDM",sred->ignore_dm,&sred->ignore_dm,0);CHKERRQ(ierr);
585   ierr = PetscOptionsTail();CHKERRQ(ierr);
586   PetscFunctionReturn(0);
587 }
588 
589 /* PC simplementation specific API's */
590 
591 static PetscErrorCode PCTelescopeGetKSP_Telescope(PC pc,KSP *ksp)
592 {
593   PC_Telescope red = (PC_Telescope)pc->data;
594   if (ksp) *ksp = red->ksp;
595   return(0);
596 }
597 
598 static PetscErrorCode PCTelescopeGetReductionFactor_Telescope(PC pc,PetscInt *fact)
599 {
600   PC_Telescope red = (PC_Telescope)pc->data;
601   if (fact) *fact = red->redfactor;
602   return(0);
603 }
604 
605 static PetscErrorCode PCTelescopeSetReductionFactor_Telescope(PC pc,PetscInt fact)
606 {
607   PC_Telescope     red = (PC_Telescope)pc->data;
608   PetscMPIInt      size;
609   PetscErrorCode   ierr;
610 
611   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
612   if (fact <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be positive",fact);
613   if (fact > size) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Reduction factor of telescoping PC %D must be <= comm.size",fact);
614   red->redfactor = fact;
615   return(0);
616 }
617 
618 static PetscErrorCode PCTelescopeGetIgnoreDM_Telescope(PC pc,PetscBool *v)
619 {
620   PC_Telescope red = (PC_Telescope)pc->data;
621   if (v) *v = red->ignore_dm;
622   return(0);
623 }
624 static PetscErrorCode PCTelescopeSetIgnoreDM_Telescope(PC pc,PetscBool v)
625 {
626   PC_Telescope red = (PC_Telescope)pc->data;
627   red->ignore_dm = v;
628   return(0);
629 }
630 
631 static PetscErrorCode PCTelescopeGetDM_Telescope(PC pc,DM *dm)
632 {
633   PC_Telescope red = (PC_Telescope)pc->data;
634   *dm = private_PCTelescopeGetSubDM(red);
635   return(0);
636 }
637 
638 /*@
639  PCTelescopeGetKSP - Gets the KSP created by the telescoping PC.
640 
641  Not Collective
642 
643  Input Parameter:
644  .  pc - the preconditioner context
645 
646  Output Parameter:
647  .  subksp - the KSP defined the smaller set of processes
648 
649  Level: advanced
650 
651  .keywords: PC, telescopting solve
652  @*/
653 PetscErrorCode PCTelescopeGetKSP(PC pc,KSP *subksp)
654 {
655   PetscTryMethod(pc,"PCTelescopeGetKSP_C",(PC,KSP*),(pc,subksp));
656   return(0);
657 }
658 
659 /*@
660  PCTelescopeGetReductionFactor - Gets the factor by which the original number of processes has been reduced by.
661 
662  Not Collective
663 
664  Input Parameter:
665  .  pc - the preconditioner context
666 
667  Output Parameter:
668  .  fact - the reduction factor
669 
670  Level: advanced
671 
672  .keywords: PC, telescoping solve
673  @*/
674 PetscErrorCode PCTelescopeGetReductionFactor(PC pc,PetscInt *fact)
675 {
676   PetscTryMethod(pc,"PCTelescopeGetReductionFactor_C",(PC,PetscInt*),(pc,fact));
677   return(0);
678 }
679 
680 /*@
681  PCTelescopeSetReductionFactor - Sets the factor by which the original number of processes has been reduced by.
682 
683  Not Collective
684 
685  Input Parameter:
686  .  pc - the preconditioner context
687 
688  Output Parameter:
689  .  fact - the reduction factor
690 
691  Level: advanced
692 
693  .keywords: PC, telescoping solve
694  @*/
695 PetscErrorCode PCTelescopeSetReductionFactor(PC pc,PetscInt fact)
696 {
697   PetscTryMethod(pc,"PCTelescopeSetReductionFactor_C",(PC,PetscInt),(pc,fact));
698   return(0);
699 }
700 
701 /*@
702  PCTelescopeGetIgnoreDM - Get the flag indicating if any DM attached to the PC will be used.
703 
704  Not Collective
705 
706  Input Parameter:
707  .  pc - the preconditioner context
708 
709  Output Parameter:
710  .  v - the flag
711 
712  Level: advanced
713 
714  .keywords: PC, telescoping solve
715  @*/
716 PetscErrorCode PCTelescopeGetIgnoreDM(PC pc,PetscBool *v)
717 {
718   PetscTryMethod(pc,"PCTelescopeGetIgnoreDM_C",(PC,PetscBool*),(pc,v));
719   return(0);
720 }
721 
722 /*@
723  PCTelescopeSetIgnoreDM - Set a flag to ignore any DM attached to the PC.
724 
725  Not Collective
726 
727  Input Parameter:
728  .  pc - the preconditioner context
729 
730  Output Parameter:
731  .  v - Use PETSC_TRUE to ignore any DM
732 
733  Level: advanced
734 
735  .keywords: PC, telescoping solve
736  @*/
737 PetscErrorCode PCTelescopeSetIgnoreDM(PC pc,PetscBool v)
738 {
739   PetscTryMethod(pc,"PCTelescopeSetIgnoreDM_C",(PC,PetscBool),(pc,v));
740   return(0);
741 }
742 
743 /*@
744  PCTelescopeGetDM - Get the re-partitioned DM attached to the sub KSP.
745 
746  Not Collective
747 
748  Input Parameter:
749  .  pc - the preconditioner context
750 
751  Output Parameter:
752  .  subdm - The re-partitioned DM
753 
754  Level: advanced
755 
756  .keywords: PC, telescoping solve
757  @*/
758 PetscErrorCode PCTelescopeGetDM(PC pc,DM *subdm)
759 {
760   PetscTryMethod(pc,"PCTelescopeGetDM_C",(PC,DM*),(pc,subdm));
761   return(0);
762 }
763 
764 /* -------------------------------------------------------------------------------------*/
765 /*MC
766    PCTELESCOPE - Runs a KSP solver on a sub-group of processors. MPI processes not in the sub-communicator are idle during the solve.
767 
768    Options Database:
769 +  -pc_telescope_reduction_factor <n> - factor to use communicator size by, for example if you are using 64 MPI processes and
770    use an n of 4, the new sub-communicator will be 4 defined with 64/4 processes
771 -  -pc_telescope_ignore_dm <false> - flag to indicate whether an attached DM should be ignored
772 
773    Level: advanced
774 
775    Notes:
776    The default KSP is PREONLY. If a DM is attached to the PC, it is re-partitioned on the sub-communicator.
777    Both the B mat operator and the right hand side vector are permuted into the new DOF ordering defined by the re-partitioned DM.
778    Currently only support for re-partitioning a DMDA is provided.
779    Any nullspace attached to the original Bmat operator are extracted, re-partitioned and set on the repartitioned Bmat operator.
780    KSPSetComputeOperators() is not propagated to the sub KSP.
781    Currently there is no support for the flag -pc_use_amat
782 
783    Optimizations:
784    (i) VecPlaceArray() could be used for scatters between the vectors defined on different sized communicators, thereby slightly reducing memory footprint;
785    (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).
786 
787  Contributed by Dave May
788 
789 .seealso:  PCTelescopeGetKSP(), PCTelescopeGetDM(),
790  PCTelescopeGetReductionFactor(), PCTelescopeSetReductionFactor(),
791  PCTelescopeGetIgnoreDM(), PCTelescopeSetIgnoreDM(), PCREDUNDANT
792 M*/
793 #undef __FUNCT__
794 #define __FUNCT__ "PCCreate_Telescope"
795 PETSC_EXTERN PetscErrorCode PCCreate_Telescope(PC pc)
796 {
797   PetscErrorCode       ierr;
798   struct _PC_Telescope *sred;
799 
800   PetscFunctionBegin;
801   ierr = PetscNewLog(pc,&sred);CHKERRQ(ierr);
802   sred->redfactor      = 1;
803   sred->ignore_dm      = PETSC_FALSE;
804   pc->data             = (void*)sred;
805 
806   pc->ops->apply          = PCApply_Telescope;
807   pc->ops->applytranspose = NULL;
808   pc->ops->setup          = PCSetUp_Telescope;
809   pc->ops->destroy        = PCDestroy_Telescope;
810   pc->ops->reset          = PCReset_Telescope;
811   pc->ops->setfromoptions = PCSetFromOptions_Telescope;
812   pc->ops->view           = PCView_Telescope;
813 
814   sred->pctelescope_setup_type              = PCTelescopeSetUp_default;
815   sred->pctelescope_matcreate_type          = PCTelescopeMatCreate_default;
816   sred->pctelescope_matnullspacecreate_type = PCTelescopeMatNullSpaceCreate_default;
817   sred->pctelescope_reset_type              = NULL;
818 
819   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetKSP_C",PCTelescopeGetKSP_Telescope);CHKERRQ(ierr);
820   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetReductionFactor_C",PCTelescopeGetReductionFactor_Telescope);CHKERRQ(ierr);
821   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetReductionFactor_C",PCTelescopeSetReductionFactor_Telescope);CHKERRQ(ierr);
822   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetIgnoreDM_C",PCTelescopeGetIgnoreDM_Telescope);CHKERRQ(ierr);
823   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeSetIgnoreDM_C",PCTelescopeSetIgnoreDM_Telescope);CHKERRQ(ierr);
824   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCTelescopeGetDM_C",PCTelescopeGetDM_Telescope);CHKERRQ(ierr);
825   PetscFunctionReturn(0);
826 }
827