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