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