xref: /petsc/src/ksp/pc/impls/telescope/telescope_dmda.c (revision 2e65eb737b5b07432530db55f6b2a145ebc548b2)
1 
2 #include <petsc/private/matimpl.h>
3 #include <petsc/private/pcimpl.h>
4 #include <petsc/private/dmimpl.h>
5 #include <petscksp.h>           /*I "petscksp.h" I*/
6 #include <petscdm.h>
7 #include <petscdmda.h>
8 
9 #include "../src/ksp/pc/impls/telescope/telescope.h"
10 
11 static PetscBool  cited = PETSC_FALSE;
12 static const char citation[] =
13 "@inproceedings{MaySananRuppKnepleySmith2016,\n"
14 "  title     = {Extreme-Scale Multigrid Components within PETSc},\n"
15 "  author    = {Dave A. May and Patrick Sanan and Karl Rupp and Matthew G. Knepley and Barry F. Smith},\n"
16 "  booktitle = {Proceedings of the Platform for Advanced Scientific Computing Conference},\n"
17 "  series    = {PASC '16},\n"
18 "  isbn      = {978-1-4503-4126-4},\n"
19 "  location  = {Lausanne, Switzerland},\n"
20 "  pages     = {5:1--5:12},\n"
21 "  articleno = {5},\n"
22 "  numpages  = {12},\n"
23 "  url       = {https://doi.acm.org/10.1145/2929908.2929913},\n"
24 "  doi       = {10.1145/2929908.2929913},\n"
25 "  acmid     = {2929913},\n"
26 "  publisher = {ACM},\n"
27 "  address   = {New York, NY, USA},\n"
28 "  keywords  = {GPU, HPC, agglomeration, coarse-level solver, multigrid, parallel computing, preconditioning},\n"
29 "  year      = {2016}\n"
30 "}\n";
31 
32 static PetscErrorCode _DMDADetermineRankFromGlobalIJK(PetscInt dim,PetscInt i,PetscInt j,PetscInt k,
33                                                       PetscInt Mp,PetscInt Np,PetscInt Pp,
34                                                       PetscInt start_i[],PetscInt start_j[],PetscInt start_k[],
35                                                       PetscInt span_i[],PetscInt span_j[],PetscInt span_k[],
36                                                       PetscMPIInt *_pi,PetscMPIInt *_pj,PetscMPIInt *_pk,PetscMPIInt *rank_re)
37 {
38   PetscInt pi,pj,pk,n;
39 
40   PetscFunctionBegin;
41   *rank_re = -1;
42   if (_pi) *_pi = -1;
43   if (_pj) *_pj = -1;
44   if (_pk) *_pk = -1;
45   pi = pj = pk = -1;
46   if (_pi) {
47     for (n=0; n<Mp; n++) {
48       if ((i >= start_i[n]) && (i < start_i[n]+span_i[n])) {
49         pi = n;
50         break;
51       }
52     }
53     PetscCheck(pi != -1,PETSC_COMM_SELF,PETSC_ERR_USER,"[dmda-ijk] pi cannot be determined : range %" PetscInt_FMT ", val %" PetscInt_FMT,Mp,i);
54     *_pi = pi;
55   }
56 
57   if (_pj) {
58     for (n=0; n<Np; n++) {
59       if ((j >= start_j[n]) && (j < start_j[n]+span_j[n])) {
60         pj = n;
61         break;
62       }
63     }
64     PetscCheck(pj != -1,PETSC_COMM_SELF,PETSC_ERR_USER,"[dmda-ijk] pj cannot be determined : range %" PetscInt_FMT ", val %" PetscInt_FMT,Np,j);
65     *_pj = pj;
66   }
67 
68   if (_pk) {
69     for (n=0; n<Pp; n++) {
70       if ((k >= start_k[n]) && (k < start_k[n]+span_k[n])) {
71         pk = n;
72         break;
73       }
74     }
75     PetscCheck(pk != -1,PETSC_COMM_SELF,PETSC_ERR_USER,"[dmda-ijk] pk cannot be determined : range %" PetscInt_FMT ", val %" PetscInt_FMT,Pp,k);
76     *_pk = pk;
77   }
78 
79   switch (dim) {
80   case 1:
81     *rank_re = pi;
82     break;
83   case 2:
84     *rank_re = pi + pj * Mp;
85     break;
86   case 3:
87     *rank_re = pi + pj * Mp + pk * (Mp*Np);
88     break;
89   }
90   PetscFunctionReturn(0);
91 }
92 
93 static PetscErrorCode _DMDADetermineGlobalS0(PetscInt dim,PetscMPIInt rank_re,PetscInt Mp_re,PetscInt Np_re,PetscInt Pp_re,
94                                       PetscInt range_i_re[],PetscInt range_j_re[],PetscInt range_k_re[],PetscInt *s0)
95 {
96   PetscInt i,j,k,start_IJK = 0;
97   PetscInt rank_ijk;
98 
99   PetscFunctionBegin;
100   switch (dim) {
101   case 1:
102     for (i=0; i<Mp_re; i++) {
103       rank_ijk = i;
104       if (rank_ijk < rank_re) {
105         start_IJK += range_i_re[i];
106       }
107     }
108     break;
109     case 2:
110     for (j=0; j<Np_re; j++) {
111       for (i=0; i<Mp_re; i++) {
112         rank_ijk = i + j*Mp_re;
113         if (rank_ijk < rank_re) {
114           start_IJK += range_i_re[i]*range_j_re[j];
115         }
116       }
117     }
118     break;
119     case 3:
120     for (k=0; k<Pp_re; k++) {
121       for (j=0; j<Np_re; j++) {
122         for (i=0; i<Mp_re; i++) {
123           rank_ijk = i + j*Mp_re + k*Mp_re*Np_re;
124           if (rank_ijk < rank_re) {
125             start_IJK += range_i_re[i]*range_j_re[j]*range_k_re[k];
126           }
127         }
128       }
129     }
130     break;
131   }
132   *s0 = start_IJK;
133   PetscFunctionReturn(0);
134 }
135 
136 static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors2d(PC_Telescope sred,DM dm,DM subdm)
137 {
138   DM             cdm;
139   Vec            coor,coor_natural,perm_coors;
140   PetscInt       i,j,si,sj,ni,nj,M,N,Ml,Nl,c,nidx;
141   PetscInt       *fine_indices;
142   IS             is_fine,is_local;
143   VecScatter     sctx;
144 
145   PetscFunctionBegin;
146   PetscCall(DMGetCoordinates(dm,&coor));
147   if (!coor) return(0);
148   if (PCTelescope_isActiveRank(sred)) {
149     PetscCall(DMDASetUniformCoordinates(subdm,0.0,1.0,0.0,1.0,0.0,1.0));
150   }
151   /* Get the coordinate vector from the distributed array */
152   PetscCall(DMGetCoordinateDM(dm,&cdm));
153   PetscCall(DMDACreateNaturalVector(cdm,&coor_natural));
154 
155   PetscCall(DMDAGlobalToNaturalBegin(cdm,coor,INSERT_VALUES,coor_natural));
156   PetscCall(DMDAGlobalToNaturalEnd(cdm,coor,INSERT_VALUES,coor_natural));
157 
158   /* get indices of the guys I want to grab */
159   PetscCall(DMDAGetInfo(dm,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
160   if (PCTelescope_isActiveRank(sred)) {
161     PetscCall(DMDAGetCorners(subdm,&si,&sj,NULL,&ni,&nj,NULL));
162     Ml = ni;
163     Nl = nj;
164   } else {
165     si = sj = 0;
166     ni = nj = 0;
167     Ml = Nl = 0;
168   }
169 
170   PetscCall(PetscMalloc1(Ml*Nl*2,&fine_indices));
171   c = 0;
172   if (PCTelescope_isActiveRank(sred)) {
173     for (j=sj; j<sj+nj; j++) {
174       for (i=si; i<si+ni; i++) {
175         nidx = (i) + (j)*M;
176         fine_indices[c  ] = 2 * nidx     ;
177         fine_indices[c+1] = 2 * nidx + 1 ;
178         c = c + 2;
179       }
180     }
181     PetscCheck(c == Ml*Nl*2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"c %" PetscInt_FMT " should equal 2 * Ml %" PetscInt_FMT " * Nl %" PetscInt_FMT,c,Ml,Nl);
182   }
183 
184   /* generate scatter */
185   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm),Ml*Nl*2,fine_indices,PETSC_USE_POINTER,&is_fine));
186   PetscCall(ISCreateStride(PETSC_COMM_SELF,Ml*Nl*2,0,1,&is_local));
187 
188   /* scatter */
189   PetscCall(VecCreate(PETSC_COMM_SELF,&perm_coors));
190   PetscCall(VecSetSizes(perm_coors,PETSC_DECIDE,Ml*Nl*2));
191   PetscCall(VecSetType(perm_coors,VECSEQ));
192 
193   PetscCall(VecScatterCreate(coor_natural,is_fine,perm_coors,is_local,&sctx));
194   PetscCall(VecScatterBegin(sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD));
195   PetscCall(VecScatterEnd(  sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD));
196   /* access */
197   if (PCTelescope_isActiveRank(sred)) {
198     Vec               _coors;
199     const PetscScalar *LA_perm;
200     PetscScalar       *LA_coors;
201 
202     PetscCall(DMGetCoordinates(subdm,&_coors));
203     PetscCall(VecGetArrayRead(perm_coors,&LA_perm));
204     PetscCall(VecGetArray(_coors,&LA_coors));
205     for (i=0; i<Ml*Nl*2; i++) {
206       LA_coors[i] = LA_perm[i];
207     }
208     PetscCall(VecRestoreArray(_coors,&LA_coors));
209     PetscCall(VecRestoreArrayRead(perm_coors,&LA_perm));
210   }
211 
212   /* update local coords */
213   if (PCTelescope_isActiveRank(sred)) {
214     DM  _dmc;
215     Vec _coors,_coors_local;
216     PetscCall(DMGetCoordinateDM(subdm,&_dmc));
217     PetscCall(DMGetCoordinates(subdm,&_coors));
218     PetscCall(DMGetCoordinatesLocal(subdm,&_coors_local));
219     PetscCall(DMGlobalToLocalBegin(_dmc,_coors,INSERT_VALUES,_coors_local));
220     PetscCall(DMGlobalToLocalEnd(_dmc,_coors,INSERT_VALUES,_coors_local));
221   }
222   PetscCall(VecScatterDestroy(&sctx));
223   PetscCall(ISDestroy(&is_fine));
224   PetscCall(PetscFree(fine_indices));
225   PetscCall(ISDestroy(&is_local));
226   PetscCall(VecDestroy(&perm_coors));
227   PetscCall(VecDestroy(&coor_natural));
228   PetscFunctionReturn(0);
229 }
230 
231 static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors3d(PC_Telescope sred,DM dm,DM subdm)
232 {
233   DM             cdm;
234   Vec            coor,coor_natural,perm_coors;
235   PetscInt       i,j,k,si,sj,sk,ni,nj,nk,M,N,P,Ml,Nl,Pl,c,nidx;
236   PetscInt       *fine_indices;
237   IS             is_fine,is_local;
238   VecScatter     sctx;
239 
240   PetscFunctionBegin;
241   PetscCall(DMGetCoordinates(dm,&coor));
242   if (!coor) PetscFunctionReturn(0);
243 
244   if (PCTelescope_isActiveRank(sred)) {
245     PetscCall(DMDASetUniformCoordinates(subdm,0.0,1.0,0.0,1.0,0.0,1.0));
246   }
247 
248   /* Get the coordinate vector from the distributed array */
249   PetscCall(DMGetCoordinateDM(dm,&cdm));
250   PetscCall(DMDACreateNaturalVector(cdm,&coor_natural));
251   PetscCall(DMDAGlobalToNaturalBegin(cdm,coor,INSERT_VALUES,coor_natural));
252   PetscCall(DMDAGlobalToNaturalEnd(cdm,coor,INSERT_VALUES,coor_natural));
253 
254   /* get indices of the guys I want to grab */
255   PetscCall(DMDAGetInfo(dm,NULL,&M,&N,&P,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
256 
257   if (PCTelescope_isActiveRank(sred)) {
258     PetscCall(DMDAGetCorners(subdm,&si,&sj,&sk,&ni,&nj,&nk));
259     Ml = ni;
260     Nl = nj;
261     Pl = nk;
262   } else {
263     si = sj = sk = 0;
264     ni = nj = nk = 0;
265     Ml = Nl = Pl = 0;
266   }
267 
268   PetscCall(PetscMalloc1(Ml*Nl*Pl*3,&fine_indices));
269 
270   c = 0;
271   if (PCTelescope_isActiveRank(sred)) {
272     for (k=sk; k<sk+nk; k++) {
273       for (j=sj; j<sj+nj; j++) {
274         for (i=si; i<si+ni; i++) {
275           nidx = (i) + (j)*M + (k)*M*N;
276           fine_indices[c  ] = 3 * nidx     ;
277           fine_indices[c+1] = 3 * nidx + 1 ;
278           fine_indices[c+2] = 3 * nidx + 2 ;
279           c = c + 3;
280         }
281       }
282     }
283   }
284 
285   /* generate scatter */
286   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm),Ml*Nl*Pl*3,fine_indices,PETSC_USE_POINTER,&is_fine));
287   PetscCall(ISCreateStride(PETSC_COMM_SELF,Ml*Nl*Pl*3,0,1,&is_local));
288 
289   /* scatter */
290   PetscCall(VecCreate(PETSC_COMM_SELF,&perm_coors));
291   PetscCall(VecSetSizes(perm_coors,PETSC_DECIDE,Ml*Nl*Pl*3));
292   PetscCall(VecSetType(perm_coors,VECSEQ));
293   PetscCall(VecScatterCreate(coor_natural,is_fine,perm_coors,is_local,&sctx));
294   PetscCall(VecScatterBegin(sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD));
295   PetscCall(VecScatterEnd(  sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD));
296 
297   /* access */
298   if (PCTelescope_isActiveRank(sred)) {
299     Vec               _coors;
300     const PetscScalar *LA_perm;
301     PetscScalar       *LA_coors;
302 
303     PetscCall(DMGetCoordinates(subdm,&_coors));
304     PetscCall(VecGetArrayRead(perm_coors,&LA_perm));
305     PetscCall(VecGetArray(_coors,&LA_coors));
306     for (i=0; i<Ml*Nl*Pl*3; i++) {
307       LA_coors[i] = LA_perm[i];
308     }
309     PetscCall(VecRestoreArray(_coors,&LA_coors));
310     PetscCall(VecRestoreArrayRead(perm_coors,&LA_perm));
311   }
312 
313   /* update local coords */
314   if (PCTelescope_isActiveRank(sred)) {
315     DM  _dmc;
316     Vec _coors,_coors_local;
317 
318     PetscCall(DMGetCoordinateDM(subdm,&_dmc));
319     PetscCall(DMGetCoordinates(subdm,&_coors));
320     PetscCall(DMGetCoordinatesLocal(subdm,&_coors_local));
321     PetscCall(DMGlobalToLocalBegin(_dmc,_coors,INSERT_VALUES,_coors_local));
322     PetscCall(DMGlobalToLocalEnd(_dmc,_coors,INSERT_VALUES,_coors_local));
323   }
324 
325   PetscCall(VecScatterDestroy(&sctx));
326   PetscCall(ISDestroy(&is_fine));
327   PetscCall(PetscFree(fine_indices));
328   PetscCall(ISDestroy(&is_local));
329   PetscCall(VecDestroy(&perm_coors));
330   PetscCall(VecDestroy(&coor_natural));
331   PetscFunctionReturn(0);
332 }
333 
334 static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx)
335 {
336   PetscInt       dim;
337   DM             dm,subdm;
338   PetscSubcomm   psubcomm;
339   MPI_Comm       comm;
340   Vec            coor;
341 
342   PetscFunctionBegin;
343   PetscCall(PCGetDM(pc,&dm));
344   PetscCall(DMGetCoordinates(dm,&coor));
345   if (!coor) PetscFunctionReturn(0);
346   psubcomm = sred->psubcomm;
347   comm = PetscSubcommParent(psubcomm);
348   subdm = ctx->dmrepart;
349 
350   PetscCall(PetscInfo(pc,"PCTelescope: setting up the coordinates (DMDA)\n"));
351   PetscCall(DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
352   switch (dim) {
353   case 1:
354     SETERRQ(comm,PETSC_ERR_SUP,"Telescope: DMDA (1D) repartitioning not provided");
355   case 2: PetscCall(PCTelescopeSetUp_dmda_repart_coors2d(sred,dm,subdm));
356     break;
357   case 3: PetscCall(PCTelescopeSetUp_dmda_repart_coors3d(sred,dm,subdm));
358     break;
359   }
360   PetscFunctionReturn(0);
361 }
362 
363 /* setup repartitioned dm */
364 PetscErrorCode PCTelescopeSetUp_dmda_repart(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx)
365 {
366   DM                    dm;
367   PetscInt              dim,nx,ny,nz,ndof,nsw,sum,k;
368   DMBoundaryType        bx,by,bz;
369   DMDAStencilType       stencil;
370   const PetscInt        *_range_i_re;
371   const PetscInt        *_range_j_re;
372   const PetscInt        *_range_k_re;
373   DMDAInterpolationType itype;
374   PetscInt              refine_x,refine_y,refine_z;
375   MPI_Comm              comm,subcomm;
376   const char            *prefix;
377 
378   PetscFunctionBegin;
379   comm = PetscSubcommParent(sred->psubcomm);
380   subcomm = PetscSubcommChild(sred->psubcomm);
381   PetscCall(PCGetDM(pc,&dm));
382 
383   PetscCall(DMDAGetInfo(dm,&dim,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,&nsw,&bx,&by,&bz,&stencil));
384   PetscCall(DMDAGetInterpolationType(dm,&itype));
385   PetscCall(DMDAGetRefinementFactor(dm,&refine_x,&refine_y,&refine_z));
386 
387   ctx->dmrepart = NULL;
388   _range_i_re = _range_j_re = _range_k_re = NULL;
389   /* Create DMDA on the child communicator */
390   if (PCTelescope_isActiveRank(sred)) {
391     switch (dim) {
392     case 1:
393       PetscCall(PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (1D)\n"));
394       /* PetscCall(DMDACreate1d(subcomm,bx,nx,ndof,nsw,NULL,&ctx->dmrepart)); */
395       ny = nz = 1;
396       by = bz = DM_BOUNDARY_NONE;
397       break;
398     case 2:
399       PetscCall(PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (2D)\n"));
400       /* PetscCall(DMDACreate2d(subcomm,bx,by,stencil,nx,ny, PETSC_DECIDE,PETSC_DECIDE,
401          ndof,nsw, NULL,NULL,&ctx->dmrepart)); */
402       nz = 1;
403       bz = DM_BOUNDARY_NONE;
404       break;
405     case 3:
406       PetscCall(PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (3D)\n"));
407       /* PetscCall(DMDACreate3d(subcomm,bx,by,bz,stencil,nx,ny,nz,
408          PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE, ndof,nsw, NULL,NULL,NULL,&ctx->dmrepart)); */
409       break;
410     }
411     /*
412      The API DMDACreate1d(), DMDACreate2d(), DMDACreate3d() does not allow us to set/append
413      a unique option prefix for the DM, thus I prefer to expose the contents of these API's here.
414      This allows users to control the partitioning of the subDM.
415     */
416     PetscCall(DMDACreate(subcomm,&ctx->dmrepart));
417     /* Set unique option prefix name */
418     PetscCall(KSPGetOptionsPrefix(sred->ksp,&prefix));
419     PetscCall(DMSetOptionsPrefix(ctx->dmrepart,prefix));
420     PetscCall(DMAppendOptionsPrefix(ctx->dmrepart,"repart_"));
421     /* standard setup from DMDACreate{1,2,3}d() */
422     PetscCall(DMSetDimension(ctx->dmrepart,dim));
423     PetscCall(DMDASetSizes(ctx->dmrepart,nx,ny,nz));
424     PetscCall(DMDASetNumProcs(ctx->dmrepart,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE));
425     PetscCall(DMDASetBoundaryType(ctx->dmrepart,bx,by,bz));
426     PetscCall(DMDASetDof(ctx->dmrepart,ndof));
427     PetscCall(DMDASetStencilType(ctx->dmrepart,stencil));
428     PetscCall(DMDASetStencilWidth(ctx->dmrepart,nsw));
429     PetscCall(DMDASetOwnershipRanges(ctx->dmrepart,NULL,NULL,NULL));
430     PetscCall(DMSetFromOptions(ctx->dmrepart));
431     PetscCall(DMSetUp(ctx->dmrepart));
432     /* Set refinement factors and interpolation type from the partent */
433     PetscCall(DMDASetRefinementFactor(ctx->dmrepart,refine_x,refine_y,refine_z));
434     PetscCall(DMDASetInterpolationType(ctx->dmrepart,itype));
435 
436     PetscCall(DMDAGetInfo(ctx->dmrepart,NULL,NULL,NULL,NULL,&ctx->Mp_re,&ctx->Np_re,&ctx->Pp_re,NULL,NULL,NULL,NULL,NULL,NULL));
437     PetscCall(DMDAGetOwnershipRanges(ctx->dmrepart,&_range_i_re,&_range_j_re,&_range_k_re));
438 
439     ctx->dmrepart->ops->creatematrix = dm->ops->creatematrix;
440     ctx->dmrepart->ops->createdomaindecomposition = dm->ops->createdomaindecomposition;
441   }
442 
443   /* generate ranges for repartitioned dm */
444   /* note - assume rank 0 always participates */
445   /* TODO: use a single MPI call */
446   PetscCallMPI(MPI_Bcast(&ctx->Mp_re,1,MPIU_INT,0,comm));
447   PetscCallMPI(MPI_Bcast(&ctx->Np_re,1,MPIU_INT,0,comm));
448   PetscCallMPI(MPI_Bcast(&ctx->Pp_re,1,MPIU_INT,0,comm));
449 
450   PetscCall(PetscCalloc3(ctx->Mp_re,&ctx->range_i_re,ctx->Np_re,&ctx->range_j_re,ctx->Pp_re,&ctx->range_k_re));
451 
452   if (_range_i_re) PetscCall(PetscArraycpy(ctx->range_i_re,_range_i_re,ctx->Mp_re));
453   if (_range_j_re) PetscCall(PetscArraycpy(ctx->range_j_re,_range_j_re,ctx->Np_re));
454   if (_range_k_re) PetscCall(PetscArraycpy(ctx->range_k_re,_range_k_re,ctx->Pp_re));
455 
456   /* TODO: use a single MPI call */
457   PetscCallMPI(MPI_Bcast(ctx->range_i_re,ctx->Mp_re,MPIU_INT,0,comm));
458   PetscCallMPI(MPI_Bcast(ctx->range_j_re,ctx->Np_re,MPIU_INT,0,comm));
459   PetscCallMPI(MPI_Bcast(ctx->range_k_re,ctx->Pp_re,MPIU_INT,0,comm));
460 
461   PetscCall(PetscMalloc3(ctx->Mp_re,&ctx->start_i_re,ctx->Np_re,&ctx->start_j_re,ctx->Pp_re,&ctx->start_k_re));
462 
463   sum = 0;
464   for (k=0; k<ctx->Mp_re; k++) {
465     ctx->start_i_re[k] = sum;
466     sum += ctx->range_i_re[k];
467   }
468 
469   sum = 0;
470   for (k=0; k<ctx->Np_re; k++) {
471     ctx->start_j_re[k] = sum;
472     sum += ctx->range_j_re[k];
473   }
474 
475   sum = 0;
476   for (k=0; k<ctx->Pp_re; k++) {
477     ctx->start_k_re[k] = sum;
478     sum += ctx->range_k_re[k];
479   }
480 
481   /* attach repartitioned dm to child ksp */
482   {
483     PetscErrorCode (*dmksp_func)(KSP,Mat,Mat,void*);
484     void           *dmksp_ctx;
485 
486     PetscCall(DMKSPGetComputeOperators(dm,&dmksp_func,&dmksp_ctx));
487 
488     /* attach dm to ksp on sub communicator */
489     if (PCTelescope_isActiveRank(sred)) {
490       PetscCall(KSPSetDM(sred->ksp,ctx->dmrepart));
491 
492       if (!dmksp_func || sred->ignore_kspcomputeoperators) {
493         PetscCall(KSPSetDMActive(sred->ksp,PETSC_FALSE));
494       } else {
495         /* sub ksp inherits dmksp_func and context provided by user */
496         PetscCall(KSPSetComputeOperators(sred->ksp,dmksp_func,dmksp_ctx));
497         PetscCall(KSPSetDMActive(sred->ksp,PETSC_TRUE));
498       }
499     }
500   }
501   PetscFunctionReturn(0);
502 }
503 
504 PetscErrorCode PCTelescopeSetUp_dmda_permutation_3d(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx)
505 {
506   DM             dm;
507   MPI_Comm       comm;
508   Mat            Pscalar,P;
509   PetscInt       ndof;
510   PetscInt       i,j,k,location,startI[3],endI[3],lenI[3],nx,ny,nz;
511   PetscInt       sr,er,Mr;
512   Vec            V;
513 
514   PetscFunctionBegin;
515   PetscCall(PetscInfo(pc,"PCTelescope: setting up the permutation matrix (DMDA-3D)\n"));
516   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
517 
518   PetscCall(PCGetDM(pc,&dm));
519   PetscCall(DMDAGetInfo(dm,NULL,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,NULL,NULL,NULL,NULL,NULL));
520 
521   PetscCall(DMGetGlobalVector(dm,&V));
522   PetscCall(VecGetSize(V,&Mr));
523   PetscCall(VecGetOwnershipRange(V,&sr,&er));
524   PetscCall(DMRestoreGlobalVector(dm,&V));
525   sr = sr / ndof;
526   er = er / ndof;
527   Mr = Mr / ndof;
528 
529   PetscCall(MatCreate(comm,&Pscalar));
530   PetscCall(MatSetSizes(Pscalar,(er-sr),(er-sr),Mr,Mr));
531   PetscCall(MatSetType(Pscalar,MATAIJ));
532   PetscCall(MatSeqAIJSetPreallocation(Pscalar,1,NULL));
533   PetscCall(MatMPIAIJSetPreallocation(Pscalar,1,NULL,1,NULL));
534 
535   PetscCall(DMDAGetCorners(dm,NULL,NULL,NULL,&lenI[0],&lenI[1],&lenI[2]));
536   PetscCall(DMDAGetCorners(dm,&startI[0],&startI[1],&startI[2],&endI[0],&endI[1],&endI[2]));
537   endI[0] += startI[0];
538   endI[1] += startI[1];
539   endI[2] += startI[2];
540 
541   for (k=startI[2]; k<endI[2]; k++) {
542     for (j=startI[1]; j<endI[1]; j++) {
543       for (i=startI[0]; i<endI[0]; i++) {
544         PetscMPIInt rank_ijk_re,rank_reI[3];
545         PetscInt    s0_re;
546         PetscInt    ii,jj,kk,local_ijk_re,mapped_ijk;
547         PetscInt    lenI_re[3];
548 
549         location = (i - startI[0]) + (j - startI[1])*lenI[0] + (k - startI[2])*lenI[0]*lenI[1];
550         PetscCall(_DMDADetermineRankFromGlobalIJK(3,i,j,k,   ctx->Mp_re,ctx->Np_re,ctx->Pp_re,
551                                                   ctx->start_i_re,ctx->start_j_re,ctx->start_k_re,
552                                                   ctx->range_i_re,ctx->range_j_re,ctx->range_k_re,
553                                                   &rank_reI[0],&rank_reI[1],&rank_reI[2],&rank_ijk_re));
554         PetscCall(_DMDADetermineGlobalS0(3,rank_ijk_re, ctx->Mp_re,ctx->Np_re,ctx->Pp_re, ctx->range_i_re,ctx->range_j_re,ctx->range_k_re, &s0_re));
555         ii = i - ctx->start_i_re[ rank_reI[0] ];
556         PetscCheck(ii >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error ii");
557         jj = j - ctx->start_j_re[ rank_reI[1] ];
558         PetscCheck(jj >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error jj");
559         kk = k - ctx->start_k_re[ rank_reI[2] ];
560         PetscCheck(kk >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error kk");
561         lenI_re[0] = ctx->range_i_re[ rank_reI[0] ];
562         lenI_re[1] = ctx->range_j_re[ rank_reI[1] ];
563         lenI_re[2] = ctx->range_k_re[ rank_reI[2] ];
564         local_ijk_re = ii + jj * lenI_re[0] + kk * lenI_re[0] * lenI_re[1];
565         mapped_ijk = s0_re + local_ijk_re;
566         PetscCall(MatSetValue(Pscalar,sr+location,mapped_ijk,1.0,INSERT_VALUES));
567       }
568     }
569   }
570   PetscCall(MatAssemblyBegin(Pscalar,MAT_FINAL_ASSEMBLY));
571   PetscCall(MatAssemblyEnd(Pscalar,MAT_FINAL_ASSEMBLY));
572   PetscCall(MatCreateMAIJ(Pscalar,ndof,&P));
573   PetscCall(MatDestroy(&Pscalar));
574   ctx->permutation = P;
575   PetscFunctionReturn(0);
576 }
577 
578 PetscErrorCode PCTelescopeSetUp_dmda_permutation_2d(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx)
579 {
580   DM             dm;
581   MPI_Comm       comm;
582   Mat            Pscalar,P;
583   PetscInt       ndof;
584   PetscInt       i,j,location,startI[2],endI[2],lenI[2],nx,ny,nz;
585   PetscInt       sr,er,Mr;
586   Vec            V;
587 
588   PetscFunctionBegin;
589   PetscCall(PetscInfo(pc,"PCTelescope: setting up the permutation matrix (DMDA-2D)\n"));
590   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
591   PetscCall(PCGetDM(pc,&dm));
592   PetscCall(DMDAGetInfo(dm,NULL,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,NULL,NULL,NULL,NULL,NULL));
593   PetscCall(DMGetGlobalVector(dm,&V));
594   PetscCall(VecGetSize(V,&Mr));
595   PetscCall(VecGetOwnershipRange(V,&sr,&er));
596   PetscCall(DMRestoreGlobalVector(dm,&V));
597   sr = sr / ndof;
598   er = er / ndof;
599   Mr = Mr / ndof;
600 
601   PetscCall(MatCreate(comm,&Pscalar));
602   PetscCall(MatSetSizes(Pscalar,(er-sr),(er-sr),Mr,Mr));
603   PetscCall(MatSetType(Pscalar,MATAIJ));
604   PetscCall(MatSeqAIJSetPreallocation(Pscalar,1,NULL));
605   PetscCall(MatMPIAIJSetPreallocation(Pscalar,1,NULL,1,NULL));
606 
607   PetscCall(DMDAGetCorners(dm,NULL,NULL,NULL,&lenI[0],&lenI[1],NULL));
608   PetscCall(DMDAGetCorners(dm,&startI[0],&startI[1],NULL,&endI[0],&endI[1],NULL));
609   endI[0] += startI[0];
610   endI[1] += startI[1];
611 
612   for (j=startI[1]; j<endI[1]; j++) {
613     for (i=startI[0]; i<endI[0]; i++) {
614       PetscMPIInt rank_ijk_re,rank_reI[3];
615       PetscInt    s0_re;
616       PetscInt    ii,jj,local_ijk_re,mapped_ijk;
617       PetscInt    lenI_re[3];
618 
619       location = (i - startI[0]) + (j - startI[1])*lenI[0];
620       PetscCall(_DMDADetermineRankFromGlobalIJK(2,i,j,0,   ctx->Mp_re,ctx->Np_re,ctx->Pp_re,
621                                                 ctx->start_i_re,ctx->start_j_re,ctx->start_k_re,
622                                                 ctx->range_i_re,ctx->range_j_re,ctx->range_k_re,
623                                                 &rank_reI[0],&rank_reI[1],NULL,&rank_ijk_re));
624 
625       PetscCall(_DMDADetermineGlobalS0(2,rank_ijk_re, ctx->Mp_re,ctx->Np_re,ctx->Pp_re, ctx->range_i_re,ctx->range_j_re,ctx->range_k_re, &s0_re));
626 
627       ii = i - ctx->start_i_re[ rank_reI[0] ];
628       PetscCheck(ii >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm2d] index error ii");
629       jj = j - ctx->start_j_re[ rank_reI[1] ];
630       PetscCheck(jj >= 0,PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm2d] index error jj");
631 
632       lenI_re[0] = ctx->range_i_re[ rank_reI[0] ];
633       lenI_re[1] = ctx->range_j_re[ rank_reI[1] ];
634       local_ijk_re = ii + jj * lenI_re[0];
635       mapped_ijk = s0_re + local_ijk_re;
636       PetscCall(MatSetValue(Pscalar,sr+location,mapped_ijk,1.0,INSERT_VALUES));
637     }
638   }
639   PetscCall(MatAssemblyBegin(Pscalar,MAT_FINAL_ASSEMBLY));
640   PetscCall(MatAssemblyEnd(Pscalar,MAT_FINAL_ASSEMBLY));
641   PetscCall(MatCreateMAIJ(Pscalar,ndof,&P));
642   PetscCall(MatDestroy(&Pscalar));
643   ctx->permutation = P;
644   PetscFunctionReturn(0);
645 }
646 
647 PetscErrorCode PCTelescopeSetUp_dmda_scatters(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx)
648 {
649   Vec            xred,yred,xtmp,x,xp;
650   VecScatter     scatter;
651   IS             isin;
652   Mat            B;
653   PetscInt       m,bs,st,ed;
654   MPI_Comm       comm;
655 
656   PetscFunctionBegin;
657   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
658   PetscCall(PCGetOperators(pc,NULL,&B));
659   PetscCall(MatCreateVecs(B,&x,NULL));
660   PetscCall(MatGetBlockSize(B,&bs));
661   PetscCall(VecDuplicate(x,&xp));
662   m = 0;
663   xred = NULL;
664   yred = NULL;
665   if (PCTelescope_isActiveRank(sred)) {
666     PetscCall(DMCreateGlobalVector(ctx->dmrepart,&xred));
667     PetscCall(VecDuplicate(xred,&yred));
668     PetscCall(VecGetOwnershipRange(xred,&st,&ed));
669     PetscCall(ISCreateStride(comm,ed-st,st,1,&isin));
670     PetscCall(VecGetLocalSize(xred,&m));
671   } else {
672     PetscCall(VecGetOwnershipRange(x,&st,&ed));
673     PetscCall(ISCreateStride(comm,0,st,1,&isin));
674   }
675   PetscCall(ISSetBlockSize(isin,bs));
676   PetscCall(VecCreate(comm,&xtmp));
677   PetscCall(VecSetSizes(xtmp,m,PETSC_DECIDE));
678   PetscCall(VecSetBlockSize(xtmp,bs));
679   PetscCall(VecSetType(xtmp,((PetscObject)x)->type_name));
680   PetscCall(VecScatterCreate(x,isin,xtmp,NULL,&scatter));
681   sred->xred    = xred;
682   sred->yred    = yred;
683   sred->isin    = isin;
684   sred->scatter = scatter;
685   sred->xtmp    = xtmp;
686 
687   ctx->xp       = xp;
688   PetscCall(VecDestroy(&x));
689   PetscFunctionReturn(0);
690 }
691 
692 PetscErrorCode PCTelescopeSetUp_dmda(PC pc,PC_Telescope sred)
693 {
694   PC_Telescope_DMDACtx *ctx;
695   PetscInt             dim;
696   DM                   dm;
697   MPI_Comm             comm;
698 
699   PetscFunctionBegin;
700   PetscCall(PetscInfo(pc,"PCTelescope: setup (DMDA)\n"));
701   PetscCall(PetscNew(&ctx));
702   sred->dm_ctx = (void*)ctx;
703 
704   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
705   PetscCall(PCGetDM(pc,&dm));
706   PetscCall(DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
707 
708   PCTelescopeSetUp_dmda_repart(pc,sred,ctx);
709   PCTelescopeSetUp_dmda_repart_coors(pc,sred,ctx);
710   switch (dim) {
711   case 1:
712     SETERRQ(comm,PETSC_ERR_SUP,"Telescope: DMDA (1D) repartitioning not provided");
713   case 2: PetscCall(PCTelescopeSetUp_dmda_permutation_2d(pc,sred,ctx));
714     break;
715   case 3: PetscCall(PCTelescopeSetUp_dmda_permutation_3d(pc,sred,ctx));
716     break;
717   }
718   PetscCall(PCTelescopeSetUp_dmda_scatters(pc,sred,ctx));
719   PetscFunctionReturn(0);
720 }
721 
722 PetscErrorCode PCTelescopeMatCreate_dmda_dmactivefalse(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A)
723 {
724   PC_Telescope_DMDACtx *ctx;
725   MPI_Comm             comm,subcomm;
726   Mat                  Bperm,Bred,B,P;
727   PetscInt             nr,nc;
728   IS                   isrow,iscol;
729   Mat                  Blocal,*_Blocal;
730 
731   PetscFunctionBegin;
732   PetscCall(PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (DMDA)\n"));
733   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
734   subcomm = PetscSubcommChild(sred->psubcomm);
735   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;
736 
737   PetscCall(PCGetOperators(pc,NULL,&B));
738   PetscCall(MatGetSize(B,&nr,&nc));
739 
740   P = ctx->permutation;
741   PetscCall(MatPtAP(B,P,MAT_INITIAL_MATRIX,1.1,&Bperm));
742 
743   /* Get submatrices */
744   isrow = sred->isin;
745   PetscCall(ISCreateStride(comm,nc,0,1,&iscol));
746 
747   PetscCall(MatCreateSubMatrices(Bperm,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal));
748   Blocal = *_Blocal;
749   Bred = NULL;
750   if (PCTelescope_isActiveRank(sred)) {
751     PetscInt mm;
752 
753     if (reuse != MAT_INITIAL_MATRIX) {Bred = *A;}
754     PetscCall(MatGetSize(Blocal,&mm,NULL));
755     /* PetscCall(MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,PETSC_DECIDE,reuse,&Bred)); */
756     PetscCall(MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred));
757   }
758   *A = Bred;
759 
760   PetscCall(ISDestroy(&iscol));
761   PetscCall(MatDestroy(&Bperm));
762   PetscCall(MatDestroyMatrices(1,&_Blocal));
763   PetscFunctionReturn(0);
764 }
765 
766 PetscErrorCode PCTelescopeMatCreate_dmda(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A)
767 {
768   DM             dm;
769   PetscErrorCode (*dmksp_func)(KSP,Mat,Mat,void*);
770   void           *dmksp_ctx;
771 
772   PetscFunctionBegin;
773   PetscCall(PCGetDM(pc,&dm));
774   PetscCall(DMKSPGetComputeOperators(dm,&dmksp_func,&dmksp_ctx));
775   /* We assume that dmksp_func = NULL, is equivalent to dmActive = PETSC_FALSE */
776   if (dmksp_func && !sred->ignore_kspcomputeoperators) {
777     DM  dmrepart;
778     Mat Ak;
779 
780     *A = NULL;
781     if (PCTelescope_isActiveRank(sred)) {
782       PetscCall(KSPGetDM(sred->ksp,&dmrepart));
783       if (reuse == MAT_INITIAL_MATRIX) {
784         PetscCall(DMCreateMatrix(dmrepart,&Ak));
785         *A = Ak;
786       } else if (reuse == MAT_REUSE_MATRIX) {
787         Ak = *A;
788       }
789       /*
790        There is no need to explicitly assemble the operator now,
791        the sub-KSP will call the method provided to KSPSetComputeOperators() during KSPSetUp()
792       */
793     }
794   } else {
795     PetscCall(PCTelescopeMatCreate_dmda_dmactivefalse(pc,sred,reuse,A));
796   }
797   PetscFunctionReturn(0);
798 }
799 
800 PetscErrorCode PCTelescopeSubNullSpaceCreate_dmda_Telescope(PC pc,PC_Telescope sred,MatNullSpace nullspace,MatNullSpace *sub_nullspace)
801 {
802   PetscBool            has_const;
803   PetscInt             i,k,n = 0;
804   const Vec            *vecs;
805   Vec                  *sub_vecs = NULL;
806   MPI_Comm             subcomm;
807   PC_Telescope_DMDACtx *ctx;
808 
809   PetscFunctionBegin;
810   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;
811   subcomm = PetscSubcommChild(sred->psubcomm);
812   PetscCall(MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs));
813 
814   if (PCTelescope_isActiveRank(sred)) {
815     /* create new vectors */
816     if (n) {
817       PetscCall(VecDuplicateVecs(sred->xred,n,&sub_vecs));
818     }
819   }
820 
821   /* copy entries */
822   for (k=0; k<n; k++) {
823     const PetscScalar *x_array;
824     PetscScalar       *LA_sub_vec;
825     PetscInt          st,ed;
826 
827     /* permute vector into ordering associated with re-partitioned dmda */
828     PetscCall(MatMultTranspose(ctx->permutation,vecs[k],ctx->xp));
829 
830     /* pull in vector x->xtmp */
831     PetscCall(VecScatterBegin(sred->scatter,ctx->xp,sred->xtmp,INSERT_VALUES,SCATTER_FORWARD));
832     PetscCall(VecScatterEnd(sred->scatter,ctx->xp,sred->xtmp,INSERT_VALUES,SCATTER_FORWARD));
833 
834     /* copy vector entries into xred */
835     PetscCall(VecGetArrayRead(sred->xtmp,&x_array));
836     if (sub_vecs) {
837       if (sub_vecs[k]) {
838         PetscCall(VecGetOwnershipRange(sub_vecs[k],&st,&ed));
839         PetscCall(VecGetArray(sub_vecs[k],&LA_sub_vec));
840         for (i=0; i<ed-st; i++) {
841           LA_sub_vec[i] = x_array[i];
842         }
843         PetscCall(VecRestoreArray(sub_vecs[k],&LA_sub_vec));
844       }
845     }
846     PetscCall(VecRestoreArrayRead(sred->xtmp,&x_array));
847   }
848 
849   if (PCTelescope_isActiveRank(sred)) {
850     /* create new (near) nullspace for redundant object */
851     PetscCall(MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,sub_nullspace));
852     PetscCall(VecDestroyVecs(n,&sub_vecs));
853     PetscCheck(!nullspace->remove,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope");
854     PetscCheck(!nullspace->rmctx,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callback context not supported when propagating (near) nullspaces with PCTelescope");
855   }
856   PetscFunctionReturn(0);
857 }
858 
859 PetscErrorCode PCTelescopeMatNullSpaceCreate_dmda(PC pc,PC_Telescope sred,Mat sub_mat)
860 {
861   Mat            B;
862 
863   PetscFunctionBegin;
864   PetscCall(PCGetOperators(pc,NULL,&B));
865   {
866     MatNullSpace nullspace,sub_nullspace;
867     PetscCall(MatGetNullSpace(B,&nullspace));
868     if (nullspace) {
869       PetscCall(PetscInfo(pc,"PCTelescope: generating nullspace (DMDA)\n"));
870       PetscCall(PCTelescopeSubNullSpaceCreate_dmda_Telescope(pc,sred,nullspace,&sub_nullspace));
871       if (PCTelescope_isActiveRank(sred)) {
872         PetscCall(MatSetNullSpace(sub_mat,sub_nullspace));
873         PetscCall(MatNullSpaceDestroy(&sub_nullspace));
874       }
875     }
876   }
877   {
878     MatNullSpace nearnullspace,sub_nearnullspace;
879     PetscCall(MatGetNearNullSpace(B,&nearnullspace));
880     if (nearnullspace) {
881       PetscCall(PetscInfo(pc,"PCTelescope: generating near nullspace (DMDA)\n"));
882       PetscCall(PCTelescopeSubNullSpaceCreate_dmda_Telescope(pc,sred,nearnullspace,&sub_nearnullspace));
883       if (PCTelescope_isActiveRank(sred)) {
884         PetscCall(MatSetNearNullSpace(sub_mat,sub_nearnullspace));
885         PetscCall(MatNullSpaceDestroy(&sub_nearnullspace));
886       }
887     }
888   }
889   PetscFunctionReturn(0);
890 }
891 
892 PetscErrorCode PCApply_Telescope_dmda(PC pc,Vec x,Vec y)
893 {
894   PC_Telescope         sred = (PC_Telescope)pc->data;
895   Mat                  perm;
896   Vec                  xtmp,xp,xred,yred;
897   PetscInt             i,st,ed;
898   VecScatter           scatter;
899   PetscScalar          *array;
900   const PetscScalar    *x_array;
901   PC_Telescope_DMDACtx *ctx;
902 
903   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;
904   xtmp    = sred->xtmp;
905   scatter = sred->scatter;
906   xred    = sred->xred;
907   yred    = sred->yred;
908   perm    = ctx->permutation;
909   xp      = ctx->xp;
910 
911   PetscFunctionBegin;
912   PetscCall(PetscCitationsRegister(citation,&cited));
913 
914   /* permute vector into ordering associated with re-partitioned dmda */
915   PetscCall(MatMultTranspose(perm,x,xp));
916 
917   /* pull in vector x->xtmp */
918   PetscCall(VecScatterBegin(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD));
919   PetscCall(VecScatterEnd(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD));
920 
921   /* copy vector entries into xred */
922   PetscCall(VecGetArrayRead(xtmp,&x_array));
923   if (xred) {
924     PetscScalar *LA_xred;
925     PetscCall(VecGetOwnershipRange(xred,&st,&ed));
926 
927     PetscCall(VecGetArray(xred,&LA_xred));
928     for (i=0; i<ed-st; i++) {
929       LA_xred[i] = x_array[i];
930     }
931     PetscCall(VecRestoreArray(xred,&LA_xred));
932   }
933   PetscCall(VecRestoreArrayRead(xtmp,&x_array));
934 
935   /* solve */
936   if (PCTelescope_isActiveRank(sred)) {
937     PetscCall(KSPSolve(sred->ksp,xred,yred));
938     PetscCall(KSPCheckSolve(sred->ksp,pc,yred));
939   }
940 
941   /* return vector */
942   PetscCall(VecGetArray(xtmp,&array));
943   if (yred) {
944     const PetscScalar *LA_yred;
945     PetscCall(VecGetOwnershipRange(yred,&st,&ed));
946     PetscCall(VecGetArrayRead(yred,&LA_yred));
947     for (i=0; i<ed-st; i++) {
948       array[i] = LA_yred[i];
949     }
950     PetscCall(VecRestoreArrayRead(yred,&LA_yred));
951   }
952   PetscCall(VecRestoreArray(xtmp,&array));
953   PetscCall(VecScatterBegin(scatter,xtmp,xp,INSERT_VALUES,SCATTER_REVERSE));
954   PetscCall(VecScatterEnd(scatter,xtmp,xp,INSERT_VALUES,SCATTER_REVERSE));
955   PetscCall(MatMult(perm,xp,y));
956   PetscFunctionReturn(0);
957 }
958 
959 PetscErrorCode PCApplyRichardson_Telescope_dmda(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason)
960 {
961   PC_Telescope         sred = (PC_Telescope)pc->data;
962   Mat                  perm;
963   Vec                  xtmp,xp,yred;
964   PetscInt             i,st,ed;
965   VecScatter           scatter;
966   const PetscScalar    *x_array;
967   PetscBool            default_init_guess_value = PETSC_FALSE;
968   PC_Telescope_DMDACtx *ctx;
969 
970   PetscFunctionBegin;
971   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;
972   xtmp    = sred->xtmp;
973   scatter = sred->scatter;
974   yred    = sred->yred;
975   perm    = ctx->permutation;
976   xp      = ctx->xp;
977 
978   PetscCheck(its <= 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCApplyRichardson_Telescope_dmda only supports max_it = 1");
979   *reason = (PCRichardsonConvergedReason)0;
980 
981   if (!zeroguess) {
982     PetscCall(PetscInfo(pc,"PCTelescopeDMDA: Scattering y for non-zero-initial guess\n"));
983     /* permute vector into ordering associated with re-partitioned dmda */
984     PetscCall(MatMultTranspose(perm,y,xp));
985 
986     /* pull in vector x->xtmp */
987     PetscCall(VecScatterBegin(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD));
988     PetscCall(VecScatterEnd(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD));
989 
990     /* copy vector entries into xred */
991     PetscCall(VecGetArrayRead(xtmp,&x_array));
992     if (yred) {
993       PetscScalar *LA_yred;
994       PetscCall(VecGetOwnershipRange(yred,&st,&ed));
995       PetscCall(VecGetArray(yred,&LA_yred));
996       for (i=0; i<ed-st; i++) {
997         LA_yred[i] = x_array[i];
998       }
999       PetscCall(VecRestoreArray(yred,&LA_yred));
1000     }
1001     PetscCall(VecRestoreArrayRead(xtmp,&x_array));
1002   }
1003 
1004   if (PCTelescope_isActiveRank(sred)) {
1005     PetscCall(KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value));
1006     if (!zeroguess) PetscCall(KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE));
1007   }
1008 
1009   PetscCall(PCApply_Telescope_dmda(pc,x,y));
1010 
1011   if (PCTelescope_isActiveRank(sred)) {
1012     PetscCall(KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value));
1013   }
1014 
1015   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1016   *outits = 1;
1017   PetscFunctionReturn(0);
1018 }
1019 
1020 PetscErrorCode PCReset_Telescope_dmda(PC pc)
1021 {
1022   PC_Telescope         sred = (PC_Telescope)pc->data;
1023   PC_Telescope_DMDACtx *ctx;
1024 
1025   PetscFunctionBegin;
1026   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;
1027   PetscCall(VecDestroy(&ctx->xp));
1028   PetscCall(MatDestroy(&ctx->permutation));
1029   PetscCall(DMDestroy(&ctx->dmrepart));
1030   PetscCall(PetscFree3(ctx->range_i_re,ctx->range_j_re,ctx->range_k_re));
1031   PetscCall(PetscFree3(ctx->start_i_re,ctx->start_j_re,ctx->start_k_re));
1032   PetscFunctionReturn(0);
1033 }
1034 
1035 PetscErrorCode DMView_DA_Short_3d(DM dm,PetscViewer v)
1036 {
1037   PetscInt       M,N,P,m,n,p,ndof,nsw;
1038   MPI_Comm       comm;
1039   PetscMPIInt    size;
1040   const char*    prefix;
1041 
1042   PetscFunctionBegin;
1043   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
1044   PetscCallMPI(MPI_Comm_size(comm,&size));
1045   PetscCall(DMGetOptionsPrefix(dm,&prefix));
1046   PetscCall(DMDAGetInfo(dm,NULL,&M,&N,&P,&m,&n,&p,&ndof,&nsw,NULL,NULL,NULL,NULL));
1047   if (prefix) PetscCall(PetscViewerASCIIPrintf(v,"DMDA Object:    (%s)    %d MPI processes\n",prefix,size));
1048   else PetscCall(PetscViewerASCIIPrintf(v,"DMDA Object:    %d MPI processes\n",size));
1049   PetscCall(PetscViewerASCIIPrintf(v,"  M %" PetscInt_FMT " N %" PetscInt_FMT " P %" PetscInt_FMT " m %" PetscInt_FMT " n %" PetscInt_FMT " p %" PetscInt_FMT " dof %" PetscInt_FMT " overlap %" PetscInt_FMT "\n",M,N,P,m,n,p,ndof,nsw));
1050   PetscFunctionReturn(0);
1051 }
1052 
1053 PetscErrorCode DMView_DA_Short_2d(DM dm,PetscViewer v)
1054 {
1055   PetscInt       M,N,m,n,ndof,nsw;
1056   MPI_Comm       comm;
1057   PetscMPIInt    size;
1058   const char*    prefix;
1059 
1060   PetscFunctionBegin;
1061   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
1062   PetscCallMPI(MPI_Comm_size(comm,&size));
1063   PetscCall(DMGetOptionsPrefix(dm,&prefix));
1064   PetscCall(DMDAGetInfo(dm,NULL,&M,&N,NULL,&m,&n,NULL,&ndof,&nsw,NULL,NULL,NULL,NULL));
1065   if (prefix) PetscCall(PetscViewerASCIIPrintf(v,"DMDA Object:    (%s)    %d MPI processes\n",prefix,size));
1066   else PetscCall(PetscViewerASCIIPrintf(v,"DMDA Object:    %d MPI processes\n",size));
1067   PetscCall(PetscViewerASCIIPrintf(v,"  M %" PetscInt_FMT " N %" PetscInt_FMT " m %" PetscInt_FMT " n %" PetscInt_FMT " dof %" PetscInt_FMT " overlap %" PetscInt_FMT "\n",M,N,m,n,ndof,nsw));
1068   PetscFunctionReturn(0);
1069 }
1070 
1071 PetscErrorCode DMView_DA_Short(DM dm,PetscViewer v)
1072 {
1073   PetscInt       dim;
1074 
1075   PetscFunctionBegin;
1076   PetscCall(DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
1077   switch (dim) {
1078   case 2: PetscCall(DMView_DA_Short_2d(dm,v));
1079     break;
1080   case 3: PetscCall(DMView_DA_Short_3d(dm,v));
1081     break;
1082   }
1083   PetscFunctionReturn(0);
1084 }
1085