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