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