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