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