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