Actual source code: telescope_dmda.c

petsc-3.9.4 2018-09-11
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       = {http://doi.acm.org/10.1145/2929908.2929913},\n"
 25: "  doi       = {10.1145/2929908.2929913},\n"
 26: "  acmid     = {2929913},\n"
 27: "  publisher = {ACM},\n"
 28: "  address   = {New York, NY, USA},\n"
 29: "  keywords  = {GPU, HPC, agglomeration, coarse-level solver, multigrid, parallel computing, preconditioning},\n"
 30: "  year      = {2016}\n"
 31: "}\n";

 33: PetscErrorCode _DMDADetermineRankFromGlobalIJK(PetscInt dim,PetscInt i,PetscInt j,PetscInt k,
 34:                                                PetscInt Mp,PetscInt Np,PetscInt Pp,
 35:                                                PetscInt start_i[],PetscInt start_j[],PetscInt start_k[],
 36:                                                PetscInt span_i[],PetscInt span_j[],PetscInt span_k[],
 37:                                                PetscMPIInt *_pi,PetscMPIInt *_pj,PetscMPIInt *_pk,PetscMPIInt *rank_re)
 38: {
 39:   PetscInt pi,pj,pk,n;

 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: 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: PetscErrorCode PCTelescopeSetUp_dmda_repart_coors2d(PetscSubcomm psubcomm,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 (isActiveRank(psubcomm)) {
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 (isActiveRank(psubcomm)) {
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 (isActiveRank(psubcomm)) {
175:     for (j=sj; j<sj+nj; j++) {
176:       for (i=si; i<si+ni; i++) {
177:         nidx = (i) + (j)*M;
178:         fine_indices[c  ] = 2 * nidx     ;
179:         fine_indices[c+1] = 2 * nidx + 1 ;
180:         c = c + 2;
181:       }
182:     }
183:     if (c != Ml*Nl*2) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c %D should equal 2 * Ml %D * Nl %D",c,Ml,Nl);
184:   }

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 (isActiveRank(psubcomm)) {
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 (isActiveRank(psubcomm)) {
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: PetscErrorCode PCTelescopeSetUp_dmda_repart_coors3d(PetscSubcomm psubcomm,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 (isActiveRank(psubcomm)) {
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 (isActiveRank(psubcomm)) {
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 (isActiveRank(psubcomm)) {
275:     for (k=sk; k<sk+nk; k++) {
276:       for (j=sj; j<sj+nj; j++) {
277:         for (i=si; i<si+ni; i++) {
278:           nidx = (i) + (j)*M + (k)*M*N;
279:           fine_indices[c  ] = 3 * nidx     ;
280:           fine_indices[c+1] = 3 * nidx + 1 ;
281:           fine_indices[c+2] = 3 * nidx + 2 ;
282:           c = c + 3;
283:         }
284:       }
285:     }
286:   }

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 (isActiveRank(psubcomm)) {
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 (isActiveRank(psubcomm)) {
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: 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;


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

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

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

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

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

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

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

447:   /* generate ranges for repartitioned dm */
448:   /* note - assume rank 0 always participates */
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:   PetscCalloc1(ctx->Mp_re,&ctx->range_i_re);
454:   PetscCalloc1(ctx->Np_re,&ctx->range_j_re);
455:   PetscCalloc1(ctx->Pp_re,&ctx->range_k_re);

457:   if (_range_i_re) {PetscMemcpy(ctx->range_i_re,_range_i_re,sizeof(PetscInt)*ctx->Mp_re);}
458:   if (_range_j_re) {PetscMemcpy(ctx->range_j_re,_range_j_re,sizeof(PetscInt)*ctx->Np_re);}
459:   if (_range_k_re) {PetscMemcpy(ctx->range_k_re,_range_k_re,sizeof(PetscInt)*ctx->Pp_re);}

461:   MPI_Bcast(ctx->range_i_re,ctx->Mp_re,MPIU_INT,0,comm);
462:   MPI_Bcast(ctx->range_j_re,ctx->Np_re,MPIU_INT,0,comm);
463:   MPI_Bcast(ctx->range_k_re,ctx->Pp_re,MPIU_INT,0,comm);

465:   PetscMalloc1(ctx->Mp_re,&ctx->start_i_re);
466:   PetscMalloc1(ctx->Np_re,&ctx->start_j_re);
467:   PetscMalloc1(ctx->Pp_re,&ctx->start_k_re);

469:   sum = 0;
470:   for (k=0; k<ctx->Mp_re; k++) {
471:     ctx->start_i_re[k] = sum;
472:     sum += ctx->range_i_re[k];
473:   }

475:   sum = 0;
476:   for (k=0; k<ctx->Np_re; k++) {
477:     ctx->start_j_re[k] = sum;
478:     sum += ctx->range_j_re[k];
479:   }

481:   sum = 0;
482:   for (k=0; k<ctx->Pp_re; k++) {
483:     ctx->start_k_re[k] = sum;
484:     sum += ctx->range_k_re[k];
485:   }

487:   /* attach repartitioned dm to child ksp */
488:   {
489:     PetscErrorCode (*dmksp_func)(KSP,Mat,Mat,void*);
490:     void           *dmksp_ctx;

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

494:     /* attach dm to ksp on sub communicator */
495:     if (isActiveRank(sred->psubcomm)) {
496:       KSPSetDM(sred->ksp,ctx->dmrepart);

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

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

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

525:   PCGetDM(pc,&dm);
526:   DMDAGetInfo(dm,NULL,&nx,&ny,&nz,NULL,NULL,NULL,&ndof,NULL,NULL,NULL,NULL,NULL);

528:   DMGetGlobalVector(dm,&V);
529:   VecGetSize(V,&Mr);
530:   VecGetOwnershipRange(V,&sr,&er);
531:   DMRestoreGlobalVector(dm,&V);
532:   sr = sr / ndof;
533:   er = er / ndof;
534:   Mr = Mr / ndof;

536:   MatCreate(comm,&Pscalar);
537:   MatSetSizes(Pscalar,(er-sr),(er-sr),Mr,Mr);
538:   MatSetType(Pscalar,MATAIJ);
539:   MatSeqAIJSetPreallocation(Pscalar,1,NULL);
540:   MatMPIAIJSetPreallocation(Pscalar,1,NULL,1,NULL);

542:   DMDAGetCorners(dm,NULL,NULL,NULL,&lenI[0],&lenI[1],&lenI[2]);
543:   DMDAGetCorners(dm,&startI[0],&startI[1],&startI[2],&endI[0],&endI[1],&endI[2]);
544:   endI[0] += startI[0];
545:   endI[1] += startI[1];
546:   endI[2] += startI[2];

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

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

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

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

609:   MatCreate(comm,&Pscalar);
610:   MatSetSizes(Pscalar,(er-sr),(er-sr),Mr,Mr);
611:   MatSetType(Pscalar,MATAIJ);
612:   MatSeqAIJSetPreallocation(Pscalar,1,NULL);
613:   MatMPIAIJSetPreallocation(Pscalar,1,NULL,1,NULL);

615:   DMDAGetCorners(dm,NULL,NULL,NULL,&lenI[0],&lenI[1],NULL);
616:   DMDAGetCorners(dm,&startI[0],&startI[1],NULL,&endI[0],&endI[1],NULL);
617:   endI[0] += startI[0];
618:   endI[1] += startI[1];

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

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

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

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

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

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

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

696:   ctx->xp       = xp;
697:   VecDestroy(&x);
698:   return(0);
699: }

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

710:   PetscInfo(pc,"PCTelescope: setup (DMDA)\n");
711:   PetscMalloc1(1,&ctx);
712:   PetscMemzero(ctx,sizeof(PC_Telescope_DMDACtx));
713:   sred->dm_ctx = (void*)ctx;

715:   PetscObjectGetComm((PetscObject)pc,&comm);
716:   PCGetDM(pc,&dm);
717:   DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);

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

733: PetscErrorCode PCTelescopeMatCreate_dmda_dmactivefalse(PC pc,PC_Telescope sred,MatReuse reuse,Mat *A)
734: {
736:   PC_Telescope_DMDACtx *ctx;
737:   MPI_Comm       comm,subcomm;
738:   Mat            Bperm,Bred,B,P;
739:   PetscInt       nr,nc;
740:   IS             isrow,iscol;
741:   Mat            Blocal,*_Blocal;

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

749:   PCGetOperators(pc,NULL,&B);
750:   MatGetSize(B,&nr,&nc);

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

755:   /* Get submatrices */
756:   isrow = sred->isin;
757:   ISCreateStride(comm,nc,0,1,&iscol);

759:   MatCreateSubMatrices(Bperm,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&_Blocal);
760:   Blocal = *_Blocal;
761:   Bred = NULL;
762:   if (isActiveRank(sred->psubcomm)) {
763:     PetscInt mm;

765:     if (reuse != MAT_INITIAL_MATRIX) {Bred = *A;}
766:     MatGetSize(Blocal,&mm,NULL);
767:     /* MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,PETSC_DECIDE,reuse,&Bred); */
768:     MatCreateMPIMatConcatenateSeqMat(subcomm,Blocal,mm,reuse,&Bred);
769:   }
770:   *A = Bred;

772:   ISDestroy(&iscol);
773:   MatDestroy(&Bperm);
774:   MatDestroyMatrices(1,&_Blocal);
775:   return(0);
776: }

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

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

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

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

824:   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;
825:   subcomm = PetscSubcommChild(sred->psubcomm);
826:   MatNullSpaceGetVecs(nullspace,&has_const,&n,&vecs);

828:   if (isActiveRank(sred->psubcomm)) {
829:     /* create new vectors */
830:     if (n) {
831:       VecDuplicateVecs(sred->xred,n,&sub_vecs);
832:     }
833:   }

835:   /* copy entries */
836:   for (k=0; k<n; k++) {
837:     const PetscScalar *x_array;
838:     PetscScalar       *LA_sub_vec;
839:     PetscInt          st,ed;

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

844:     /* pull in vector x->xtmp */
845:     VecScatterBegin(sred->scatter,ctx->xp,sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);
846:     VecScatterEnd(sred->scatter,ctx->xp,sred->xtmp,INSERT_VALUES,SCATTER_FORWARD);

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

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

871:   return(0);
872: }

874: PetscErrorCode PCTelescopeMatNullSpaceCreate_dmda(PC pc,PC_Telescope sred,Mat sub_mat)
875: {
876:   PetscErrorCode   ierr;
877:   Mat              B;

880:   PCGetOperators(pc,NULL,&B);

882:   {
883:     MatNullSpace nullspace,sub_nullspace;
884:     MatGetNullSpace(B,&nullspace);
885:     if (nullspace) {
886:       PetscInfo(pc,"PCTelescope: generating nullspace (DMDA)\n");
887:       PCTelescopeSubNullSpaceCreate_dmda_Telescope(pc,sred,nullspace,&sub_nullspace);
888:       if (isActiveRank(sred->psubcomm)) {
889:         MatSetNullSpace(sub_mat,sub_nullspace);
890:         MatNullSpaceDestroy(&sub_nullspace);
891:       }
892:     }
893:   }

895:   {
896:     MatNullSpace nearnullspace,sub_nearnullspace;
897:     MatGetNullSpace(B,&nearnullspace);
898:     if (nearnullspace) {
899:       PetscInfo(pc,"PCTelescope: generating near nullspace (DMDA)\n");
900:       PCTelescopeSubNullSpaceCreate_dmda_Telescope(pc,sred,nearnullspace,&sub_nearnullspace);
901:       if (isActiveRank(sred->psubcomm)) {
902:         MatSetNullSpace(sub_mat,sub_nearnullspace);
903:         MatNullSpaceDestroy(&sub_nearnullspace);
904:       }
905:     }
906:   }
907:   return(0);
908: }

910: PetscErrorCode PCApply_Telescope_dmda(PC pc,Vec x,Vec y)
911: {
912:   PC_Telescope      sred = (PC_Telescope)pc->data;
913:   PetscErrorCode    ierr;
914:   Mat               perm;
915:   Vec               xtmp,xp,xred,yred;
916:   PetscInt          i,st,ed;
917:   VecScatter        scatter;
918:   PetscScalar       *array;
919:   const PetscScalar *x_array;
920:   PC_Telescope_DMDACtx *ctx;

922:   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;
923:   xtmp    = sred->xtmp;
924:   scatter = sred->scatter;
925:   xred    = sred->xred;
926:   yred    = sred->yred;
927:   perm  = ctx->permutation;
928:   xp    = ctx->xp;

931:   PetscCitationsRegister(citation,&cited);

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

936:   /* pull in vector x->xtmp */
937:   VecScatterBegin(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);
938:   VecScatterEnd(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);

940:   /* copy vector entires into xred */
941:   VecGetArrayRead(xtmp,&x_array);
942:   if (xred) {
943:     PetscScalar *LA_xred;
944:     VecGetOwnershipRange(xred,&st,&ed);

946:     VecGetArray(xred,&LA_xred);
947:     for (i=0; i<ed-st; i++) {
948:       LA_xred[i] = x_array[i];
949:     }
950:     VecRestoreArray(xred,&LA_xred);
951:   }
952:   VecRestoreArrayRead(xtmp,&x_array);

954:   /* solve */
955:   if (isActiveRank(sred->psubcomm)) {
956:     KSPSolve(sred->ksp,xred,yred);
957:   }

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

977: PetscErrorCode PCApplyRichardson_Telescope_dmda(PC pc,Vec x,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason)
978: {
979:   PC_Telescope      sred = (PC_Telescope)pc->data;
980:   PetscErrorCode    ierr;
981:   Mat               perm;
982:   Vec               xtmp,xp,yred;
983:   PetscInt          i,st,ed;
984:   VecScatter        scatter;
985:   const PetscScalar *x_array;
986:   PetscBool         default_init_guess_value = PETSC_FALSE;
987:   PC_Telescope_DMDACtx *ctx;

989:   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;
990:   xtmp    = sred->xtmp;
991:   scatter = sred->scatter;
992:   yred    = sred->yred;
993:   perm  = ctx->permutation;
994:   xp    = ctx->xp;

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

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

1004:     /* pull in vector x->xtmp */
1005:     VecScatterBegin(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);
1006:     VecScatterEnd(scatter,xp,xtmp,INSERT_VALUES,SCATTER_FORWARD);

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

1022:   if (isActiveRank(sred->psubcomm)) {
1023:     KSPGetInitialGuessNonzero(sred->ksp,&default_init_guess_value);
1024:     if (!zeroguess) KSPSetInitialGuessNonzero(sred->ksp,PETSC_TRUE);
1025:   }

1027:   PCApply_Telescope_dmda(pc,x,y);

1029:   if (isActiveRank(sred->psubcomm)) {
1030:     KSPSetInitialGuessNonzero(sred->ksp,default_init_guess_value);
1031:   }

1033:   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1034:   *outits = 1;
1035:   return(0);
1036: }

1038: PetscErrorCode PCReset_Telescope_dmda(PC pc)
1039: {
1040:   PetscErrorCode       ierr;
1041:   PC_Telescope         sred = (PC_Telescope)pc->data;
1042:   PC_Telescope_DMDACtx *ctx;

1045:   ctx = (PC_Telescope_DMDACtx*)sred->dm_ctx;
1046:   VecDestroy(&ctx->xp);
1047:   MatDestroy(&ctx->permutation);
1048:   DMDestroy(&ctx->dmrepart);
1049:   PetscFree(ctx->range_i_re);
1050:   PetscFree(ctx->range_j_re);
1051:   PetscFree(ctx->range_k_re);
1052:   PetscFree(ctx->start_i_re);
1053:   PetscFree(ctx->start_j_re);
1054:   PetscFree(ctx->start_k_re);
1055:   return(0);
1056: }

1058: PetscErrorCode DMView_DA_Short_3d(DM dm,PetscViewer v)
1059: {
1060:   PetscInt       M,N,P,m,n,p,ndof,nsw;
1061:   MPI_Comm       comm;
1062:   PetscMPIInt    size;
1063:   const char*    prefix;

1067:   PetscObjectGetComm((PetscObject)dm,&comm);
1068:   MPI_Comm_size(comm,&size);
1069:   DMGetOptionsPrefix(dm,&prefix);
1070:   DMDAGetInfo(dm,NULL,&M,&N,&P,&m,&n,&p,&ndof,&nsw,NULL,NULL,NULL,NULL);
1071:   if (prefix) {PetscViewerASCIIPrintf(v,"DMDA Object:    (%s)    %d MPI processes\n",prefix,size);}
1072:   else {PetscViewerASCIIPrintf(v,"DMDA Object:    %d MPI processes\n",size);}
1073:   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);
1074:   return(0);
1075: }

1077: PetscErrorCode DMView_DA_Short_2d(DM dm,PetscViewer v)
1078: {
1079:   PetscInt       M,N,m,n,ndof,nsw;
1080:   MPI_Comm       comm;
1081:   PetscMPIInt    size;
1082:   const char*    prefix;

1086:   PetscObjectGetComm((PetscObject)dm,&comm);
1087:   MPI_Comm_size(comm,&size);
1088:   DMGetOptionsPrefix(dm,&prefix);
1089:   DMDAGetInfo(dm,NULL,&M,&N,NULL,&m,&n,NULL,&ndof,&nsw,NULL,NULL,NULL,NULL);
1090:   if (prefix) {PetscViewerASCIIPrintf(v,"DMDA Object:    (%s)    %d MPI processes\n",prefix,size);}
1091:   else {PetscViewerASCIIPrintf(v,"DMDA Object:    %d MPI processes\n",size);}
1092:   PetscViewerASCIIPrintf(v,"  M %D N %D m %D n %D dof %D overlap %D\n",M,N,m,n,ndof,nsw);
1093:   return(0);
1094: }

1096: PetscErrorCode DMView_DA_Short(DM dm,PetscViewer v)
1097: {
1099:   PetscInt       dim;

1102:   DMDAGetInfo(dm,&dim,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
1103:   switch (dim) {
1104:   case 2: DMView_DA_Short_2d(dm,v);
1105:     break;
1106:   case 3: DMView_DA_Short_3d(dm,v);
1107:     break;
1108:   }
1109:   return(0);
1110: }