Actual source code: telescope_dmda.c

petsc-3.12.5 2020-03-29
Report Typos and Errors


  3:  #include <petsc/private/matimpl.h>
  4:  #include <petsc/private/pcimpl.h>
  5:  #include <petsc/private/dmimpl.h>
  6:  #include <petscksp.h>
  7:  #include <petscdm.h>
  8:  #include <petscdmda.h>

 10: #include "../src/ksp/pc/impls/telescope/telescope.h"

 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       = {https://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";

 33: static 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;

 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:   }

 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:   }

 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:   }

 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:   return(0);
 92: }

 94: static 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;

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:   return(0);
135: }

137: static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors2d(PC_Telescope sred,DM dm,DM subdm)
138: {
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;

148:   DMGetCoordinates(dm,&coor);
149:   if (!coor) return(0);
150:   if (PCTelescope_isActiveRank(sred)) {
151:     DMDASetUniformCoordinates(subdm,0.0,1.0,0.0,1.0,0.0,1.0);
152:   }
153:   /* Get the coordinate vector from the distributed array */
154:   DMGetCoordinateDM(dm,&cdm);
155:   DMDACreateNaturalVector(cdm,&coor_natural);

157:   DMDAGlobalToNaturalBegin(cdm,coor,INSERT_VALUES,coor_natural);
158:   DMDAGlobalToNaturalEnd(cdm,coor,INSERT_VALUES,coor_natural);

160:   /* get indices of the guys I want to grab */
161:   DMDAGetInfo(dm,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
162:   if (PCTelescope_isActiveRank(sred)) {
163:     DMDAGetCorners(subdm,&si,&sj,NULL,&ni,&nj,NULL);
164:     Ml = ni;
165:     Nl = nj;
166:   } else {
167:     si = sj = 0;
168:     ni = nj = 0;
169:     Ml = Nl = 0;
170:   }

172:   PetscMalloc1(Ml*Nl*2,&fine_indices);
173:   c = 0;
174:   if (PCTelescope_isActiveRank(sred)) {
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:   }

186:   /* generate scatter */
187:   ISCreateGeneral(PetscObjectComm((PetscObject)dm),Ml*Nl*2,fine_indices,PETSC_USE_POINTER,&is_fine);
188:   ISCreateStride(PETSC_COMM_SELF,Ml*Nl*2,0,1,&is_local);

190:   /* scatter */
191:   VecCreate(PETSC_COMM_SELF,&perm_coors);
192:   VecSetSizes(perm_coors,PETSC_DECIDE,Ml*Nl*2);
193:   VecSetType(perm_coors,VECSEQ);

195:   VecScatterCreate(coor_natural,is_fine,perm_coors,is_local,&sctx);
196:   VecScatterBegin(sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);
197:   VecScatterEnd(  sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);
198:   /* access */
199:   if (PCTelescope_isActiveRank(sred)) {
200:     Vec               _coors;
201:     const PetscScalar *LA_perm;
202:     PetscScalar       *LA_coors;

204:     DMGetCoordinates(subdm,&_coors);
205:     VecGetArrayRead(perm_coors,&LA_perm);
206:     VecGetArray(_coors,&LA_coors);
207:     for (i=0; i<Ml*Nl*2; i++) {
208:       LA_coors[i] = LA_perm[i];
209:     }
210:     VecRestoreArray(_coors,&LA_coors);
211:     VecRestoreArrayRead(perm_coors,&LA_perm);
212:   }

214:   /* update local coords */
215:   if (PCTelescope_isActiveRank(sred)) {
216:     DM  _dmc;
217:     Vec _coors,_coors_local;
218:     DMGetCoordinateDM(subdm,&_dmc);
219:     DMGetCoordinates(subdm,&_coors);
220:     DMGetCoordinatesLocal(subdm,&_coors_local);
221:     DMGlobalToLocalBegin(_dmc,_coors,INSERT_VALUES,_coors_local);
222:     DMGlobalToLocalEnd(_dmc,_coors,INSERT_VALUES,_coors_local);
223:   }
224:   VecScatterDestroy(&sctx);
225:   ISDestroy(&is_fine);
226:   PetscFree(fine_indices);
227:   ISDestroy(&is_local);
228:   VecDestroy(&perm_coors);
229:   VecDestroy(&coor_natural);
230:   return(0);
231: }

233: static PetscErrorCode PCTelescopeSetUp_dmda_repart_coors3d(PC_Telescope sred,DM dm,DM subdm)
234: {
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;

244:   DMGetCoordinates(dm,&coor);
245:   if (!coor) return(0);

247:   if (PCTelescope_isActiveRank(sred)) {
248:     DMDASetUniformCoordinates(subdm,0.0,1.0,0.0,1.0,0.0,1.0);
249:   }

251:   /* Get the coordinate vector from the distributed array */
252:   DMGetCoordinateDM(dm,&cdm);
253:   DMDACreateNaturalVector(cdm,&coor_natural);
254:   DMDAGlobalToNaturalBegin(cdm,coor,INSERT_VALUES,coor_natural);
255:   DMDAGlobalToNaturalEnd(cdm,coor,INSERT_VALUES,coor_natural);

257:   /* get indices of the guys I want to grab */
258:   DMDAGetInfo(dm,NULL,&M,&N,&P,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);

260:   if (PCTelescope_isActiveRank(sred)) {
261:     DMDAGetCorners(subdm,&si,&sj,&sk,&ni,&nj,&nk);
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:   }

271:   PetscMalloc1(Ml*Nl*Pl*3,&fine_indices);

273:   c = 0;
274:   if (PCTelescope_isActiveRank(sred)) {
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:   }

288:   /* generate scatter */
289:   ISCreateGeneral(PetscObjectComm((PetscObject)dm),Ml*Nl*Pl*3,fine_indices,PETSC_USE_POINTER,&is_fine);
290:   ISCreateStride(PETSC_COMM_SELF,Ml*Nl*Pl*3,0,1,&is_local);

292:   /* scatter */
293:   VecCreate(PETSC_COMM_SELF,&perm_coors);
294:   VecSetSizes(perm_coors,PETSC_DECIDE,Ml*Nl*Pl*3);
295:   VecSetType(perm_coors,VECSEQ);
296:   VecScatterCreate(coor_natural,is_fine,perm_coors,is_local,&sctx);
297:   VecScatterBegin(sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);
298:   VecScatterEnd(  sctx,coor_natural,perm_coors,INSERT_VALUES,SCATTER_FORWARD);

300:   /* access */
301:   if (PCTelescope_isActiveRank(sred)) {
302:     Vec               _coors;
303:     const PetscScalar *LA_perm;
304:     PetscScalar       *LA_coors;

306:     DMGetCoordinates(subdm,&_coors);
307:     VecGetArrayRead(perm_coors,&LA_perm);
308:     VecGetArray(_coors,&LA_coors);
309:     for (i=0; i<Ml*Nl*Pl*3; i++) {
310:       LA_coors[i] = LA_perm[i];
311:     }
312:     VecRestoreArray(_coors,&LA_coors);
313:     VecRestoreArrayRead(perm_coors,&LA_perm);
314:   }

316:   /* update local coords */
317:   if (PCTelescope_isActiveRank(sred)) {
318:     DM  _dmc;
319:     Vec _coors,_coors_local;

321:     DMGetCoordinateDM(subdm,&_dmc);
322:     DMGetCoordinates(subdm,&_coors);
323:     DMGetCoordinatesLocal(subdm,&_coors_local);
324:     DMGlobalToLocalBegin(_dmc,_coors,INSERT_VALUES,_coors_local);
325:     DMGlobalToLocalEnd(_dmc,_coors,INSERT_VALUES,_coors_local);
326:   }

328:   VecScatterDestroy(&sctx);
329:   ISDestroy(&is_fine);
330:   PetscFree(fine_indices);
331:   ISDestroy(&is_local);
332:   VecDestroy(&perm_coors);
333:   VecDestroy(&coor_natural);
334:   return(0);
335: }

337: static 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;

347:   PCGetDM(pc,&dm);
348:   DMGetCoordinates(dm,&coor);
349:   if (!coor) return(0);
350:   psubcomm = sred->psubcomm;
351:   comm = PetscSubcommParent(psubcomm);
352:   subdm = ctx->dmrepart;

354:   PetscInfo(pc,"PCTelescope: setting up the coordinates (DMDA)\n");
355:   DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
356:   switch (dim) {
357:   case 1: SETERRQ(comm,PETSC_ERR_SUP,"Telescope: DMDA (1D) repartitioning not provided");
358:     break;
359:   case 2: PCTelescopeSetUp_dmda_repart_coors2d(sred,dm,subdm);
360:     break;
361:   case 3: PCTelescopeSetUp_dmda_repart_coors3d(sred,dm,subdm);
362:     break;
363:   }
364:   return(0);
365: }

367: /* setup repartitioned dm */
368: PetscErrorCode PCTelescopeSetUp_dmda_repart(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx)
369: {
370:   PetscErrorCode        ierr;
371:   DM                    dm;
372:   PetscInt              dim,nx,ny,nz,ndof,nsw,sum,k;
373:   DMBoundaryType        bx,by,bz;
374:   DMDAStencilType       stencil;
375:   const PetscInt        *_range_i_re;
376:   const PetscInt        *_range_j_re;
377:   const PetscInt        *_range_k_re;
378:   DMDAInterpolationType itype;
379:   PetscInt              refine_x,refine_y,refine_z;
380:   MPI_Comm              comm,subcomm;
381:   const char            *prefix;

384:   comm = PetscSubcommParent(sred->psubcomm);
385:   subcomm = PetscSubcommChild(sred->psubcomm);
386:   PCGetDM(pc,&dm);

388:   DMDAGetInfo(dm,&dim,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,&nsw,&bx,&by,&bz,&stencil);
389:   DMDAGetInterpolationType(dm,&itype);
390:   DMDAGetRefinementFactor(dm,&refine_x,&refine_y,&refine_z);

392:   ctx->dmrepart = NULL;
393:   _range_i_re = _range_j_re = _range_k_re = NULL;
394:   /* Create DMDA on the child communicator */
395:   if (PCTelescope_isActiveRank(sred)) {
396:     switch (dim) {
397:     case 1:
398:       PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (1D)\n");
399:       /*DMDACreate1d(subcomm,bx,nx,ndof,nsw,NULL,&ctx->dmrepart);*/
400:       ny = nz = 1;
401:       by = bz = DM_BOUNDARY_NONE;
402:       break;
403:     case 2:
404:       PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (2D)\n");
405:       /*DMDACreate2d(subcomm,bx,by,stencil,nx,ny, PETSC_DECIDE,PETSC_DECIDE, ndof,nsw, NULL,NULL,&ctx->dmrepart);*/
406:       nz = 1;
407:       bz = DM_BOUNDARY_NONE;
408:       break;
409:     case 3:
410:       PetscInfo(pc,"PCTelescope: setting up the DMDA on comm subset (3D)\n");
411:       /*DMDACreate3d(subcomm,bx,by,bz,stencil,nx,ny,nz, PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE, ndof,nsw, NULL,NULL,NULL,&ctx->dmrepart);*/
412:       break;
413:     }
414:     /*
415:      The API DMDACreate1d(), DMDACreate2d(), DMDACreate3d() does not allow us to set/append
416:      a unique option prefix for the DM, thus I prefer to expose the contents of these API's here.
417:      This allows users to control the partitioning of the subDM.
418:     */
419:     DMDACreate(subcomm,&ctx->dmrepart);
420:     /* Set unique option prefix name */
421:     KSPGetOptionsPrefix(sred->ksp,&prefix);
422:     DMSetOptionsPrefix(ctx->dmrepart,prefix);
423:     DMAppendOptionsPrefix(ctx->dmrepart,"repart_");
424:     /* standard setup from DMDACreate{1,2,3}d() */
425:     DMSetDimension(ctx->dmrepart,dim);
426:     DMDASetSizes(ctx->dmrepart,nx,ny,nz);
427:     DMDASetNumProcs(ctx->dmrepart,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE);
428:     DMDASetBoundaryType(ctx->dmrepart,bx,by,bz);
429:     DMDASetDof(ctx->dmrepart,ndof);
430:     DMDASetStencilType(ctx->dmrepart,stencil);
431:     DMDASetStencilWidth(ctx->dmrepart,nsw);
432:     DMDASetOwnershipRanges(ctx->dmrepart,NULL,NULL,NULL);
433:     DMSetFromOptions(ctx->dmrepart);
434:     DMSetUp(ctx->dmrepart);
435:     /* Set refinement factors and interpolation type from the partent */
436:     DMDASetRefinementFactor(ctx->dmrepart,refine_x,refine_y,refine_z);
437:     DMDASetInterpolationType(ctx->dmrepart,itype);

439:     DMDAGetInfo(ctx->dmrepart,NULL,NULL,NULL,NULL,&ctx->Mp_re,&ctx->Np_re,&ctx->Pp_re,NULL,NULL,NULL,NULL,NULL,NULL);
440:     DMDAGetOwnershipRanges(ctx->dmrepart,&_range_i_re,&_range_j_re,&_range_k_re);

442:     ctx->dmrepart->ops->creatematrix = dm->ops->creatematrix;
443:     ctx->dmrepart->ops->createdomaindecomposition = dm->ops->createdomaindecomposition;
444:   }

446:   /* generate ranges for repartitioned dm */
447:   /* note - assume rank 0 always participates */
448:   /* TODO: use a single MPI call */
449:   MPI_Bcast(&ctx->Mp_re,1,MPIU_INT,0,comm);
450:   MPI_Bcast(&ctx->Np_re,1,MPIU_INT,0,comm);
451:   MPI_Bcast(&ctx->Pp_re,1,MPIU_INT,0,comm);

453:   PetscCalloc3(ctx->Mp_re,&ctx->range_i_re,ctx->Np_re,&ctx->range_j_re,ctx->Pp_re,&ctx->range_k_re);

455:   if (_range_i_re) {PetscArraycpy(ctx->range_i_re,_range_i_re,ctx->Mp_re);}
456:   if (_range_j_re) {PetscArraycpy(ctx->range_j_re,_range_j_re,ctx->Np_re);}
457:   if (_range_k_re) {PetscArraycpy(ctx->range_k_re,_range_k_re,ctx->Pp_re);}

459:   /* TODO: use a single MPI call */
460:   MPI_Bcast(ctx->range_i_re,ctx->Mp_re,MPIU_INT,0,comm);
461:   MPI_Bcast(ctx->range_j_re,ctx->Np_re,MPIU_INT,0,comm);
462:   MPI_Bcast(ctx->range_k_re,ctx->Pp_re,MPIU_INT,0,comm);

464:   PetscMalloc3(ctx->Mp_re,&ctx->start_i_re,ctx->Np_re,&ctx->start_j_re,ctx->Pp_re,&ctx->start_k_re);

466:   sum = 0;
467:   for (k=0; k<ctx->Mp_re; k++) {
468:     ctx->start_i_re[k] = sum;
469:     sum += ctx->range_i_re[k];
470:   }

472:   sum = 0;
473:   for (k=0; k<ctx->Np_re; k++) {
474:     ctx->start_j_re[k] = sum;
475:     sum += ctx->range_j_re[k];
476:   }

478:   sum = 0;
479:   for (k=0; k<ctx->Pp_re; k++) {
480:     ctx->start_k_re[k] = sum;
481:     sum += ctx->range_k_re[k];
482:   }

484:   /* attach repartitioned dm to child ksp */
485:   {
486:     PetscErrorCode (*dmksp_func)(KSP,Mat,Mat,void*);
487:     void           *dmksp_ctx;

489:     DMKSPGetComputeOperators(dm,&dmksp_func,&dmksp_ctx);

491:     /* attach dm to ksp on sub communicator */
492:     if (PCTelescope_isActiveRank(sred)) {
493:       KSPSetDM(sred->ksp,ctx->dmrepart);

495:       if (!dmksp_func || sred->ignore_kspcomputeoperators) {
496:         KSPSetDMActive(sred->ksp,PETSC_FALSE);
497:       } else {
498:         /* sub ksp inherits dmksp_func and context provided by user */
499:         KSPSetComputeOperators(sred->ksp,dmksp_func,dmksp_ctx);
500:         KSPSetDMActive(sred->ksp,PETSC_TRUE);
501:       }
502:     }
503:   }
504:   return(0);
505: }

507: PetscErrorCode PCTelescopeSetUp_dmda_permutation_3d(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx)
508: {
510:   DM             dm;
511:   MPI_Comm       comm;
512:   Mat            Pscalar,P;
513:   PetscInt       ndof;
514:   PetscInt       i,j,k,location,startI[3],endI[3],lenI[3],nx,ny,nz;
515:   PetscInt       sr,er,Mr;
516:   Vec            V;

519:   PetscInfo(pc,"PCTelescope: setting up the permutation matrix (DMDA-3D)\n");
520:   PetscObjectGetComm((PetscObject)pc,&comm);

522:   PCGetDM(pc,&dm);
523:   DMDAGetInfo(dm,NULL,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,NULL,NULL,NULL,NULL,NULL);

525:   DMGetGlobalVector(dm,&V);
526:   VecGetSize(V,&Mr);
527:   VecGetOwnershipRange(V,&sr,&er);
528:   DMRestoreGlobalVector(dm,&V);
529:   sr = sr / ndof;
530:   er = er / ndof;
531:   Mr = Mr / ndof;

533:   MatCreate(comm,&Pscalar);
534:   MatSetSizes(Pscalar,(er-sr),(er-sr),Mr,Mr);
535:   MatSetType(Pscalar,MATAIJ);
536:   MatSeqAIJSetPreallocation(Pscalar,1,NULL);
537:   MatMPIAIJSetPreallocation(Pscalar,1,NULL,1,NULL);

539:   DMDAGetCorners(dm,NULL,NULL,NULL,&lenI[0],&lenI[1],&lenI[2]);
540:   DMDAGetCorners(dm,&startI[0],&startI[1],&startI[2],&endI[0],&endI[1],&endI[2]);
541:   endI[0] += startI[0];
542:   endI[1] += startI[1];
543:   endI[2] += startI[2];

545:   for (k=startI[2]; k<endI[2]; k++) {
546:     for (j=startI[1]; j<endI[1]; j++) {
547:       for (i=startI[0]; i<endI[0]; i++) {
548:         PetscMPIInt rank_ijk_re,rank_reI[3];
549:         PetscInt    s0_re;
550:         PetscInt    ii,jj,kk,local_ijk_re,mapped_ijk;
551:         PetscInt    lenI_re[3];

553:         location = (i - startI[0]) + (j - startI[1])*lenI[0] + (k - startI[2])*lenI[0]*lenI[1];
554:         _DMDADetermineRankFromGlobalIJK(3,i,j,k,   ctx->Mp_re,ctx->Np_re,ctx->Pp_re,
555:                                                ctx->start_i_re,ctx->start_j_re,ctx->start_k_re,
556:                                                ctx->range_i_re,ctx->range_j_re,ctx->range_k_re,
557:                                                &rank_reI[0],&rank_reI[1],&rank_reI[2],&rank_ijk_re);
558:         _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);
559:         ii = i - ctx->start_i_re[ rank_reI[0] ];
560:         if (ii < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error ii");
561:         jj = j - ctx->start_j_re[ rank_reI[1] ];
562:         if (jj < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error jj");
563:         kk = k - ctx->start_k_re[ rank_reI[2] ];
564:         if (kk < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm3d] index error kk");
565:         lenI_re[0] = ctx->range_i_re[ rank_reI[0] ];
566:         lenI_re[1] = ctx->range_j_re[ rank_reI[1] ];
567:         lenI_re[2] = ctx->range_k_re[ rank_reI[2] ];
568:         local_ijk_re = ii + jj * lenI_re[0] + kk * lenI_re[0] * lenI_re[1];
569:         mapped_ijk = s0_re + local_ijk_re;
570:         MatSetValue(Pscalar,sr+location,mapped_ijk,1.0,INSERT_VALUES);
571:       }
572:     }
573:   }
574:   MatAssemblyBegin(Pscalar,MAT_FINAL_ASSEMBLY);
575:   MatAssemblyEnd(Pscalar,MAT_FINAL_ASSEMBLY);
576:   MatCreateMAIJ(Pscalar,ndof,&P);
577:   MatDestroy(&Pscalar);
578:   ctx->permutation = P;
579:   return(0);
580: }

582: PetscErrorCode PCTelescopeSetUp_dmda_permutation_2d(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx)
583: {
585:   DM             dm;
586:   MPI_Comm       comm;
587:   Mat            Pscalar,P;
588:   PetscInt       ndof;
589:   PetscInt       i,j,location,startI[2],endI[2],lenI[2],nx,ny,nz;
590:   PetscInt       sr,er,Mr;
591:   Vec            V;

594:   PetscInfo(pc,"PCTelescope: setting up the permutation matrix (DMDA-2D)\n");
595:   PetscObjectGetComm((PetscObject)pc,&comm);
596:   PCGetDM(pc,&dm);
597:   DMDAGetInfo(dm,NULL,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,NULL,NULL,NULL,NULL,NULL);
598:   DMGetGlobalVector(dm,&V);
599:   VecGetSize(V,&Mr);
600:   VecGetOwnershipRange(V,&sr,&er);
601:   DMRestoreGlobalVector(dm,&V);
602:   sr = sr / ndof;
603:   er = er / ndof;
604:   Mr = Mr / ndof;

606:   MatCreate(comm,&Pscalar);
607:   MatSetSizes(Pscalar,(er-sr),(er-sr),Mr,Mr);
608:   MatSetType(Pscalar,MATAIJ);
609:   MatSeqAIJSetPreallocation(Pscalar,1,NULL);
610:   MatMPIAIJSetPreallocation(Pscalar,1,NULL,1,NULL);

612:   DMDAGetCorners(dm,NULL,NULL,NULL,&lenI[0],&lenI[1],NULL);
613:   DMDAGetCorners(dm,&startI[0],&startI[1],NULL,&endI[0],&endI[1],NULL);
614:   endI[0] += startI[0];
615:   endI[1] += startI[1];

617:   for (j=startI[1]; j<endI[1]; j++) {
618:     for (i=startI[0]; i<endI[0]; i++) {
619:       PetscMPIInt rank_ijk_re,rank_reI[3];
620:       PetscInt    s0_re;
621:       PetscInt    ii,jj,local_ijk_re,mapped_ijk;
622:       PetscInt    lenI_re[3];

624:       location = (i - startI[0]) + (j - startI[1])*lenI[0];
625:       _DMDADetermineRankFromGlobalIJK(2,i,j,0,   ctx->Mp_re,ctx->Np_re,ctx->Pp_re,
626:                                              ctx->start_i_re,ctx->start_j_re,ctx->start_k_re,
627:                                              ctx->range_i_re,ctx->range_j_re,ctx->range_k_re,
628:                                              &rank_reI[0],&rank_reI[1],NULL,&rank_ijk_re);

630:       _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);

632:       ii = i - ctx->start_i_re[ rank_reI[0] ];
633:       if (ii < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm2d] index error ii");
634:       jj = j - ctx->start_j_re[ rank_reI[1] ];
635:       if (jj < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"[dmdarepart-perm2d] index error jj");

637:       lenI_re[0] = ctx->range_i_re[ rank_reI[0] ];
638:       lenI_re[1] = ctx->range_j_re[ rank_reI[1] ];
639:       local_ijk_re = ii + jj * lenI_re[0];
640:       mapped_ijk = s0_re + local_ijk_re;
641:       MatSetValue(Pscalar,sr+location,mapped_ijk,1.0,INSERT_VALUES);
642:     }
643:   }
644:   MatAssemblyBegin(Pscalar,MAT_FINAL_ASSEMBLY);
645:   MatAssemblyEnd(Pscalar,MAT_FINAL_ASSEMBLY);
646:   MatCreateMAIJ(Pscalar,ndof,&P);
647:   MatDestroy(&Pscalar);
648:   ctx->permutation = P;
649:   return(0);
650: }

652: PetscErrorCode PCTelescopeSetUp_dmda_scatters(PC pc,PC_Telescope sred,PC_Telescope_DMDACtx *ctx)
653: {
655:   Vec            xred,yred,xtmp,x,xp;
656:   VecScatter     scatter;
657:   IS             isin;
658:   Mat            B;
659:   PetscInt       m,bs,st,ed;
660:   MPI_Comm       comm;

663:   PetscObjectGetComm((PetscObject)pc,&comm);
664:   PCGetOperators(pc,NULL,&B);
665:   MatCreateVecs(B,&x,NULL);
666:   MatGetBlockSize(B,&bs);
667:   VecDuplicate(x,&xp);
668:   m = 0;
669:   xred = NULL;
670:   yred = NULL;
671:   if (PCTelescope_isActiveRank(sred)) {
672:     DMCreateGlobalVector(ctx->dmrepart,&xred);
673:     VecDuplicate(xred,&yred);
674:     VecGetOwnershipRange(xred,&st,&ed);
675:     ISCreateStride(comm,ed-st,st,1,&isin);
676:     VecGetLocalSize(xred,&m);
677:   } else {
678:     VecGetOwnershipRange(x,&st,&ed);
679:     ISCreateStride(comm,0,st,1,&isin);
680:   }
681:   ISSetBlockSize(isin,bs);
682:   VecCreate(comm,&xtmp);
683:   VecSetSizes(xtmp,m,PETSC_DECIDE);
684:   VecSetBlockSize(xtmp,bs);
685:   VecSetType(xtmp,((PetscObject)x)->type_name);
686:   VecScatterCreate(x,isin,xtmp,NULL,&scatter);
687:   sred->xred    = xred;
688:   sred->yred    = yred;
689:   sred->isin    = isin;
690:   sred->scatter = scatter;
691:   sred->xtmp    = xtmp;

693:   ctx->xp       = xp;
694:   VecDestroy(&x);
695:   return(0);
696: }

698: PetscErrorCode PCTelescopeSetUp_dmda(PC pc,PC_Telescope sred)
699: {
700:   PC_Telescope_DMDACtx *ctx;
701:   PetscInt             dim;
702:   DM                   dm;
703:   MPI_Comm             comm;
704:   PetscErrorCode       ierr;

707:   PetscInfo(pc,"PCTelescope: setup (DMDA)\n");
708:   PetscNew(&ctx);
709:   sred->dm_ctx = (void*)ctx;

711:   PetscObjectGetComm((PetscObject)pc,&comm);
712:   PCGetDM(pc,&dm);
713:   DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);

715:   PCTelescopeSetUp_dmda_repart(pc,sred,ctx);
716:   PCTelescopeSetUp_dmda_repart_coors(pc,sred,ctx);
717:   switch (dim) {
718:   case 1: SETERRQ(comm,PETSC_ERR_SUP,"Telescope: DMDA (1D) repartitioning not provided");
719:     break;
720:   case 2: PCTelescopeSetUp_dmda_permutation_2d(pc,sred,ctx);
721:     break;
722:   case 3: PCTelescopeSetUp_dmda_permutation_3d(pc,sred,ctx);
723:     break;
724:   }
725:   PCTelescopeSetUp_dmda_scatters(pc,sred,ctx);
726:   return(0);
727: }

729: PetscErrorCode PCTelescopeMatCreate_dmda_dmactivefalse(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A)
730: {
731:   PetscErrorCode       ierr;
732:   PC_Telescope_DMDACtx *ctx;
733:   MPI_Comm             comm,subcomm;
734:   Mat                  Bperm,Bred,B,P;
735:   PetscInt             nr,nc;
736:   IS                   isrow,iscol;
737:   Mat                  Blocal,*_Blocal;

740:   PetscInfo(pc,"PCTelescope: updating the redundant preconditioned operator (DMDA)\n");
741:   PetscObjectGetComm((PetscObject)pc,&comm);
742:   subcomm = PetscSubcommChild(sred->psubcomm);
743:   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;

745:   PCGetOperators(pc,NULL,&B);
746:   MatGetSize(B,&nr,&nc);

748:   P = ctx->permutation;
749:   MatPtAP(B,P,MAT_INITIAL_MATRIX,1.1,&Bperm);

751:   /* Get submatrices */
752:   isrow = sred->isin;
753:   ISCreateStride(comm,nc,0,1,&iscol);

755:   MatCreateSubMatrices(Bperm,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);
756:   Blocal = *_Blocal;
757:   Bred = NULL;
758:   if (PCTelescope_isActiveRank(sred)) {
759:     PetscInt mm;

761:     if (reuse != MAT_INITIAL_MATRIX) {Bred = *A;}
762:     MatGetSize(Blocal,&mm,NULL);
763:     /* MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,PETSC_DECIDE,reuse,&Bred); */
764:     MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);
765:   }
766:   *A = Bred;

768:   ISDestroy(&iscol);
769:   MatDestroy(&Bperm);
770:   MatDestroyMatrices(1,&_Blocal);
771:   return(0);
772: }

774: PetscErrorCode PCTelescopeMatCreate_dmda(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A)
775: {
777:   DM             dm;
778:   PetscErrorCode (*dmksp_func)(KSP,Mat,Mat,void*);
779:   void           *dmksp_ctx;

782:   PCGetDM(pc,&dm);
783:   DMKSPGetComputeOperators(dm,&dmksp_func,&dmksp_ctx);
784:   /* We assume that dmksp_func = NULL, is equivalent to dmActive = PETSC_FALSE */
785:   if (dmksp_func && !sred->ignore_kspcomputeoperators) {
786:     DM  dmrepart;
787:     Mat Ak;

789:     *A = NULL;
790:     if (PCTelescope_isActiveRank(sred)) {
791:       KSPGetDM(sred->ksp,&dmrepart);
792:       if (reuse == MAT_INITIAL_MATRIX) {
793:         DMCreateMatrix(dmrepart,&Ak);
794:         *A = Ak;
795:       } else if (reuse == MAT_REUSE_MATRIX) {
796:         Ak = *A;
797:       }
798:       /*
799:        There is no need to explicitly assemble the operator now,
800:        the sub-KSP will call the method provided to KSPSetComputeOperators() during KSPSetUp()
801:       */
802:     }
803:   } else {
804:     PCTelescopeMatCreate_dmda_dmactivefalse(pc,sred,reuse,A);
805:   }
806:   return(0);
807: }

809: PetscErrorCode PCTelescopeSubNullSpaceCreate_dmda_Telescope(PC pc,PC_Telescope sred,MatNullSpace nullspace,MatNullSpace *sub_nullspace)
810: {
811:   PetscErrorCode       ierr;
812:   PetscBool            has_const;
813:   PetscInt             i,k,n = 0;
814:   const Vec            *vecs;
815:   Vec                  *sub_vecs = NULL;
816:   MPI_Comm             subcomm;
817:   PC_Telescope_DMDACtx *ctx;

820:   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;
821:   subcomm = PetscSubcommChild(sred->psubcomm);
822:   MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);

824:   if (PCTelescope_isActiveRank(sred)) {
825:     /* create new vectors */
826:     if (n) {
827:       VecDuplicateVecs(sred->xred,n,&sub_vecs);
828:     }
829:   }

831:   /* copy entries */
832:   for (k=0; k<n; k++) {
833:     const PetscScalar *x_array;
834:     PetscScalar       *LA_sub_vec;
835:     PetscInt          st,ed;

837:     /* permute vector into ordering associated with re-partitioned dmda */
838:     MatMultTranspose(ctx->permutation,vecs[k],ctx->xp);

840:     /* pull in vector x->xtmp */
841:     VecScatterBegin(sred->scatter,ctx->xp,sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);
842:     VecScatterEnd(sred->scatter,ctx->xp,sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);

844:     /* copy vector entries into xred */
845:     VecGetArrayRead(sred->xtmp,&x_array);
846:     if (sub_vecs) {
847:       if (sub_vecs[k]) {
848:         VecGetOwnershipRange(sub_vecs[k],&st,&ed);
849:         VecGetArray(sub_vecs[k],&LA_sub_vec);
850:         for (i=0; i<ed-st; i++) {
851:           LA_sub_vec[i] = x_array[i];
852:         }
853:         VecRestoreArray(sub_vecs[k],&LA_sub_vec);
854:       }
855:     }
856:     VecRestoreArrayRead(sred->xtmp,&x_array);
857:   }

859:   if (PCTelescope_isActiveRank(sred)) {
860:     /* create new (near) nullspace for redundant object */
861:     MatNullSpaceCreate(subcomm,has_const,n,sub_vecs,sub_nullspace);
862:     VecDestroyVecs(n,&sub_vecs);
863:     if (nullspace->remove) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callbacks not supported when propagating (near) nullspaces with PCTelescope");
864:     if (nullspace->rmctx) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Propagation of custom remove callback context not supported when propagating (near) nullspaces with PCTelescope");
865:   }
866:   return(0);
867: }

869: PetscErrorCode PCTelescopeMatNullSpaceCreate_dmda(PC pc,PC_Telescope sred,Mat sub_mat)
870: {
872:   Mat            B;

875:   PCGetOperators(pc,NULL,&B);
876:   {
877:     MatNullSpace nullspace,sub_nullspace;
878:     MatGetNullSpace(B,&nullspace);
879:     if (nullspace) {
880:       PetscInfo(pc,"PCTelescope: generating nullspace (DMDA)\n");
881:       PCTelescopeSubNullSpaceCreate_dmda_Telescope(pc,sred,nullspace,&sub_nullspace);
882:       if (PCTelescope_isActiveRank(sred)) {
883:         MatSetNullSpace(sub_mat,sub_nullspace);
884:         MatNullSpaceDestroy(&sub_nullspace);
885:       }
886:     }
887:   }
888:   {
889:     MatNullSpace nearnullspace,sub_nearnullspace;
890:     MatGetNearNullSpace(B,&nearnullspace);
891:     if (nearnullspace) {
892:       PetscInfo(pc,"PCTelescope: generating near nullspace (DMDA)\n");
893:       PCTelescopeSubNullSpaceCreate_dmda_Telescope(pc,sred,nearnullspace,&sub_nearnullspace);
894:       if (PCTelescope_isActiveRank(sred)) {
895:         MatSetNearNullSpace(sub_mat,sub_nearnullspace);
896:         MatNullSpaceDestroy(&sub_nearnullspace);
897:       }
898:     }
899:   }
900:   return(0);
901: }

903: PetscErrorCode PCApply_Telescope_dmda(PC pc,Vec x,Vec y)
904: {
905:   PC_Telescope         sred = (PC_Telescope)pc->data;
906:   PetscErrorCode       ierr;
907:   Mat                  perm;
908:   Vec                  xtmp,xp,xred,yred;
909:   PetscInt             i,st,ed;
910:   VecScatter           scatter;
911:   PetscScalar          *array;
912:   const PetscScalar    *x_array;
913:   PC_Telescope_DMDACtx *ctx;

915:   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;
916:   xtmp    = sred->xtmp;
917:   scatter = sred->scatter;
918:   xred    = sred->xred;
919:   yred    = sred->yred;
920:   perm    = ctx->permutation;
921:   xp      = ctx->xp;

924:   PetscCitationsRegister(citation,&cited);

926:   /* permute vector into ordering associated with re-partitioned dmda */
927:   MatMultTranspose(perm,x,xp);

929:   /* pull in vector x->xtmp */
930:   VecScatterBegin(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);
931:   VecScatterEnd(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);

933:   /* copy vector entires into xred */
934:   VecGetArrayRead(xtmp,&x_array);
935:   if (xred) {
936:     PetscScalar *LA_xred;
937:     VecGetOwnershipRange(xred,&st,&ed);

939:     VecGetArray(xred,&LA_xred);
940:     for (i=0; i<ed-st; i++) {
941:       LA_xred[i] = x_array[i];
942:     }
943:     VecRestoreArray(xred,&LA_xred);
944:   }
945:   VecRestoreArrayRead(xtmp,&x_array);

947:   /* solve */
948:   if (PCTelescope_isActiveRank(sred)) {
949:     KSPSolve(sred->ksp,xred,yred);
950:     KSPCheckSolve(sred->ksp,pc,yred);
951:   }

953:   /* return vector */
954:   VecGetArray(xtmp,&array);
955:   if (yred) {
956:     const PetscScalar *LA_yred;
957:     VecGetOwnershipRange(yred,&st,&ed);
958:     VecGetArrayRead(yred,&LA_yred);
959:     for (i=0; i<ed-st; i++) {
960:       array[i] = LA_yred[i];
961:     }
962:     VecRestoreArrayRead(yred,&LA_yred);
963:   }
964:   VecRestoreArray(xtmp,&array);
965:   VecScatterBegin(scatter,xtmp,xp,INSERT_VALUES,SCATTER_REVERSE);
966:   VecScatterEnd(scatter,xtmp,xp,INSERT_VALUES,SCATTER_REVERSE);
967:   MatMult(perm,xp,y);
968:   return(0);
969: }

971: 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)
972: {
973:   PC_Telescope         sred = (PC_Telescope)pc->data;
974:   PetscErrorCode       ierr;
975:   Mat                  perm;
976:   Vec                  xtmp,xp,yred;
977:   PetscInt             i,st,ed;
978:   VecScatter           scatter;
979:   const PetscScalar    *x_array;
980:   PetscBool            default_init_guess_value = PETSC_FALSE;
981:   PC_Telescope_DMDACtx *ctx;

984:   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;
985:   xtmp    = sred->xtmp;
986:   scatter = sred->scatter;
987:   yred    = sred->yred;
988:   perm    = ctx->permutation;
989:   xp      = ctx->xp;

991:   if (its > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"PCApplyRichardson_Telescope_dmda only supports max_it = 1");
992:   *reason = (PCRichardsonConvergedReason)0;

994:   if (!zeroguess) {
995:     PetscInfo(pc,"PCTelescopeDMDA: Scattering y for non-zero-initial guess\n");
996:     /* permute vector into ordering associated with re-partitioned dmda */
997:     MatMultTranspose(perm,y,xp);

999:     /* pull in vector x->xtmp */
1000:     VecScatterBegin(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);
1001:     VecScatterEnd(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);

1003:     /* copy vector entires into xred */
1004:     VecGetArrayRead(xtmp,&x_array);
1005:     if (yred) {
1006:       PetscScalar *LA_yred;
1007:       VecGetOwnershipRange(yred,&st,&ed);
1008:       VecGetArray(yred,&LA_yred);
1009:       for (i=0; i<ed-st; i++) {
1010:         LA_yred[i] = x_array[i];
1011:       }
1012:       VecRestoreArray(yred,&LA_yred);
1013:     }
1014:     VecRestoreArrayRead(xtmp,&x_array);
1015:   }

1017:   if (PCTelescope_isActiveRank(sred)) {
1018:     KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value);
1019:     if (!zeroguess) KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE);
1020:   }

1022:   PCApply_Telescope_dmda(pc,x,y);

1024:   if (PCTelescope_isActiveRank(sred)) {
1025:     KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value);
1026:   }

1028:   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1029:   *outits = 1;
1030:   return(0);
1031: }

1033: PetscErrorCode PCReset_Telescope_dmda(PC pc)
1034: {
1035:   PetscErrorCode       ierr;
1036:   PC_Telescope         sred = (PC_Telescope)pc->data;
1037:   PC_Telescope_DMDACtx *ctx;

1040:   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;
1041:   VecDestroy(&ctx->xp);
1042:   MatDestroy(&ctx->permutation);
1043:   DMDestroy(&ctx->dmrepart);
1044:   PetscFree3(ctx->range_i_re,ctx->range_j_re,ctx->range_k_re);
1045:   PetscFree3(ctx->start_i_re,ctx->start_j_re,ctx->start_k_re);
1046:   return(0);
1047: }

1049: PetscErrorCode DMView_DA_Short_3d(DM dm,PetscViewer v)
1050: {
1051:   PetscInt       M,N,P,m,n,p,ndof,nsw;
1052:   MPI_Comm       comm;
1053:   PetscMPIInt    size;
1054:   const char*    prefix;

1058:   PetscObjectGetComm((PetscObject)dm,&comm);
1059:   MPI_Comm_size(comm,&size);
1060:   DMGetOptionsPrefix(dm,&prefix);
1061:   DMDAGetInfo(dm,NULL,&M,&N,&P,&m,&n,&p,&ndof,&nsw,NULL,NULL,NULL,NULL);
1062:   if (prefix) {PetscViewerASCIIPrintf(v,"DMDA Object:    (%s)    %d MPI processes\n",prefix,size);}
1063:   else {PetscViewerASCIIPrintf(v,"DMDA Object:    %d MPI processes\n",size);}
1064:   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);
1065:   return(0);
1066: }

1068: PetscErrorCode DMView_DA_Short_2d(DM dm,PetscViewer v)
1069: {
1070:   PetscInt       M,N,m,n,ndof,nsw;
1071:   MPI_Comm       comm;
1072:   PetscMPIInt    size;
1073:   const char*    prefix;

1077:   PetscObjectGetComm((PetscObject)dm,&comm);
1078:   MPI_Comm_size(comm,&size);
1079:   DMGetOptionsPrefix(dm,&prefix);
1080:   DMDAGetInfo(dm,NULL,&M,&N,NULL,&m,&n,NULL,&ndof,&nsw,NULL,NULL,NULL,NULL);
1081:   if (prefix) {PetscViewerASCIIPrintf(v,"DMDA Object:    (%s)    %d MPI processes\n",prefix,size);}
1082:   else {PetscViewerASCIIPrintf(v,"DMDA Object:    %d MPI processes\n",size);}
1083:   PetscViewerASCIIPrintf(v,"  M %D N %D m %D n %D dof %D overlap %D\n",M,N,m,n,ndof,nsw);
1084:   return(0);
1085: }

1087: PetscErrorCode DMView_DA_Short(DM dm,PetscViewer v)
1088: {
1090:   PetscInt       dim;

1093:   DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
1094:   switch (dim) {
1095:   case 2: DMView_DA_Short_2d(dm,v);
1096:     break;
1097:   case 3: DMView_DA_Short_3d(dm,v);
1098:     break;
1099:   }
1100:   return(0);
1101: }