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