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