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