Actual source code: forest.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1:  #include <petsc/private/dmforestimpl.h>
  2:  #include <petsc/private/dmimpl.h>
  3:  #include <petsc/private/dmlabelimpl.h>
  4:  #include <petscsf.h>

  6: PetscBool DMForestPackageInitialized = PETSC_FALSE;

  8: typedef struct _DMForestTypeLink*DMForestTypeLink;

 10: struct _DMForestTypeLink
 11: {
 12:   char             *name;
 13:   DMForestTypeLink next;
 14: };

 16: DMForestTypeLink DMForestTypeList;

 18: static PetscErrorCode DMForestPackageFinalize(void)
 19: {
 20:   DMForestTypeLink oldLink, link = DMForestTypeList;
 21:   PetscErrorCode   ierr;

 24:   while (link) {
 25:     oldLink = link;
 26:     PetscFree(oldLink->name);
 27:     link    = oldLink->next;
 28:     PetscFree(oldLink);
 29:   }
 30:   return(0);
 31: }

 33: static PetscErrorCode DMForestPackageInitialize(void)
 34: {

 38:   if (DMForestPackageInitialized) return(0);
 39:   DMForestPackageInitialized = PETSC_TRUE;

 41:   DMForestRegisterType(DMFOREST);
 42:   PetscRegisterFinalize(DMForestPackageFinalize);
 43:   return(0);
 44: }

 46: /*@C
 47:   DMForestRegisterType - Registers a DMType as a subtype of DMFOREST (so that DMIsForest() will be correct)

 49:   Not Collective

 51:   Input parameter:
 52: . name - the name of the type

 54:   Level: advanced

 56: .seealso: DMFOREST, DMIsForest()
 57: @*/
 58: PetscErrorCode DMForestRegisterType(DMType name)
 59: {
 60:   DMForestTypeLink link;
 61:   PetscErrorCode   ierr;

 64:   DMForestPackageInitialize();
 65:   PetscNew(&link);
 66:   PetscStrallocpy(name,&link->name);
 67:   link->next       = DMForestTypeList;
 68:   DMForestTypeList = link;
 69:   return(0);
 70: }

 72: /*@
 73:   DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes

 75:   Not Collective

 77:   Input parameter:
 78: . dm - the DM object

 80:   Output parameter:
 81: . isForest - whether dm is a subtype of DMFOREST

 83:   Level: intermediate

 85: .seealso: DMFOREST, DMForestRegisterType()
 86: @*/
 87: PetscErrorCode DMIsForest(DM dm, PetscBool *isForest)
 88: {
 89:   DMForestTypeLink link = DMForestTypeList;
 90:   PetscErrorCode   ierr;

 93:   while (link) {
 94:     PetscBool sameType;
 95:     PetscObjectTypeCompare((PetscObject)dm,link->name,&sameType);
 96:     if (sameType) {
 97:       *isForest = PETSC_TRUE;
 98:       return(0);
 99:     }
100:     link = link->next;
101:   }
102:   *isForest = PETSC_FALSE;
103:   return(0);
104: }

106: /*@
107:   DMForestTemplate - Create a new DM that will be adapted from a source DM.  The new DM reproduces the configuration
108:   of the source, but is not yet setup, so that the user can then define only the ways that the new DM should differ
109:   (by, e.g., refinement or repartitioning).  The source DM is also set as the adaptivity source DM of the new DM (see
110:   DMForestSetAdaptivityForest()).

112:   Collective on dm

114:   Input Parameters:
115: + dm - the source DM object
116: - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen())

118:   Output Parameter:
119: . tdm - the new DM object

121:   Level: intermediate

123: .seealso: DMForestSetAdaptivityForest()
124: @*/
125: PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm)
126: {
127:   DM_Forest                  *forest = (DM_Forest*) dm->data;
128:   DMType                     type;
129:   DM                         base;
130:   DMForestTopology           topology;
131:   PetscInt                   dim, overlap, ref, factor;
132:   DMForestAdaptivityStrategy strat;
133:   PetscDS                    ds;
134:   void                       *ctx;
135:   PetscErrorCode             (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void*);
136:   void                       *mapCtx;
137:   PetscErrorCode             ierr;

141:   DMCreate(PetscObjectComm((PetscObject)dm),tdm);
142:   DMGetType(dm,&type);
143:   DMSetType(*tdm,type);
144:   DMForestGetBaseDM(dm,&base);
145:   DMForestSetBaseDM(*tdm,base);
146:   DMForestGetTopology(dm,&topology);
147:   DMForestSetTopology(*tdm,topology);
148:   DMForestGetAdjacencyDimension(dm,&dim);
149:   DMForestSetAdjacencyDimension(*tdm,dim);
150:   DMForestGetPartitionOverlap(dm,&overlap);
151:   DMForestSetPartitionOverlap(*tdm,overlap);
152:   DMForestGetMinimumRefinement(dm,&ref);
153:   DMForestSetMinimumRefinement(*tdm,ref);
154:   DMForestGetMaximumRefinement(dm,&ref);
155:   DMForestSetMaximumRefinement(*tdm,ref);
156:   DMForestGetAdaptivityStrategy(dm,&strat);
157:   DMForestSetAdaptivityStrategy(*tdm,strat);
158:   DMForestGetGradeFactor(dm,&factor);
159:   DMForestSetGradeFactor(*tdm,factor);
160:   DMForestGetBaseCoordinateMapping(dm,&map,&mapCtx);
161:   DMForestSetBaseCoordinateMapping(*tdm,map,mapCtx);
162:   if (forest->ftemplate) {
163:     (forest->ftemplate)(dm, *tdm);
164:   }
165:   DMForestSetAdaptivityForest(*tdm,dm);
166:   DMGetDS(dm,&ds);
167:   DMSetDS(*tdm,ds);
168:   DMGetApplicationContext(dm,&ctx);
169:   DMSetApplicationContext(*tdm,&ctx);
170:   {
171:     PetscBool            isper;
172:     const PetscReal      *maxCell, *L;
173:     const DMBoundaryType *bd;

175:     DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);
176:     DMSetPeriodicity(*tdm,isper,maxCell,L,bd);
177:   }
178:   DMCopyBoundary(dm,*tdm);
179:   return(0);
180: }

182: static PetscErrorCode DMInitialize_Forest(DM dm);

184: PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm)
185: {
186:   DM_Forest      *forest = (DM_Forest*) dm->data;
187:   const char     *type;

191:   forest->refct++;
192:   (*newdm)->data = forest;
193:   PetscObjectGetType((PetscObject) dm, &type);
194:   PetscObjectChangeTypeName((PetscObject) *newdm, type);
195:   DMInitialize_Forest(*newdm);
196:   return(0);
197: }

199: static PetscErrorCode DMDestroy_Forest(DM dm)
200: {
201:   DM_Forest      *forest = (DM_Forest*) dm->data;

205:   if (--forest->refct > 0) return(0);
206:   if (forest->destroy) {forest->destroy(dm);}
207:   PetscSFDestroy(&forest->cellSF);
208:   PetscSFDestroy(&forest->preCoarseToFine);
209:   PetscSFDestroy(&forest->coarseToPreFine);
210:   DMLabelDestroy(&forest->adaptLabel);
211:   PetscFree(forest->adaptStrategy);
212:   DMDestroy(&forest->base);
213:   DMDestroy(&forest->adapt);
214:   PetscFree(forest->topology);
215:   PetscFree(forest);
216:   return(0);
217: }

219: /*@C
220:   DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase.  The topology is a string (e.g.
221:   "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest during
222:   DMSetUp().

224:   Logically collective on dm

226:   Input parameters:
227: + dm - the forest
228: - topology - the topology of the forest

230:   Level: intermediate

232: .seealso(): DMForestGetTopology(), DMForestSetBaseDM()
233: @*/
234: PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology)
235: {
236:   DM_Forest      *forest = (DM_Forest*) dm->data;

241:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup");
242:   PetscFree(forest->topology);
243:   PetscStrallocpy((const char*)topology,(char**) &forest->topology);
244:   return(0);
245: }

247: /*@C
248:   DMForestGetTopology - Get a string describing the topology of a DMForest.

250:   Not collective

252:   Input parameter:
253: . dm - the forest

255:   Output parameter:
256: . topology - the topology of the forest (e.g., 'cube', 'shell')

258:   Level: intermediate

260: .seealso: DMForestSetTopology()
261: @*/
262: PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology)
263: {
264:   DM_Forest *forest = (DM_Forest*) dm->data;

269:   *topology = forest->topology;
270:   return(0);
271: }

273: /*@
274:   DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest.  The
275:   forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its
276:   base.  In general, two forest must share a base to be comparable, to do things like construct interpolators.

278:   Logically collective on dm

280:   Input Parameters:
281: + dm - the forest
282: - base - the base DM of the forest

284:   Notes: Currently the base DM must be a DMPLEX

286:   Level: intermediate

288: .seealso(): DMForestGetBaseDM()
289: @*/
290: PetscErrorCode DMForestSetBaseDM(DM dm, DM base)
291: {
292:   DM_Forest      *forest = (DM_Forest*) dm->data;
293:   PetscInt       dim, dimEmbed;

298:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup");
299:   PetscObjectReference((PetscObject)base);
300:   DMDestroy(&forest->base);
301:   forest->base = base;
302:   if (base) {
304:     DMGetDimension(base,&dim);
305:     DMSetDimension(dm,dim);
306:     DMGetCoordinateDim(base,&dimEmbed);
307:     DMSetCoordinateDim(dm,dimEmbed);
308:   }
309:   return(0);
310: }

312: /*@
313:   DMForestGetBaseDM - Get the base DM of a DMForest forest.  The forest will be hierarchically refined from the base,
314:   and all refinements/coarsenings of the forest will share its base.  In general, two forest must share a base to be
315:   comparable, to do things like construct interpolators.

317:   Not collective

319:   Input Parameter:
320: . dm - the forest

322:   Output Parameter:
323: . base - the base DM of the forest

325:   Level: intermediate

327: .seealso(); DMForestSetBaseDM()
328: @*/
329: PetscErrorCode DMForestGetBaseDM(DM dm, DM *base)
330: {
331:   DM_Forest *forest = (DM_Forest*) dm->data;

336:   *base = forest->base;
337:   return(0);
338: }

340: PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
341: {
342:   DM_Forest *forest = (DM_Forest*) dm->data;

346:   forest->mapcoordinates    = func;
347:   forest->mapcoordinatesctx = ctx;
348:   return(0);
349: }

351: PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func) (DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
352: {
353:   DM_Forest *forest = (DM_Forest*) dm->data;

357:   if (func) *func = forest->mapcoordinates;
358:   if (ctx) *((void**) ctx) = forest->mapcoordinatesctx;
359:   return(0);
360: }

362: /*@
363:   DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be
364:   adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp().  Usually not needed
365:   by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this
366:   routine.

368:   Note that this can be called after setup with adapt = NULL, which will clear all internal data related to the
369:   adaptivity forest from dm.  This way, repeatedly adapting does not leave stale DM objects in memory.

371:   Logically collective on dm

373:   Input Parameter:
374: + dm - the new forest, which will be constructed from adapt
375: - adapt - the old forest

377:   Level: intermediate

379: .seealso: DMForestGetAdaptivityForest(), DMForestSetAdaptivityPurpose()
380: @*/
381: PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt)
382: {
383:   DM_Forest      *forest, *adaptForest, *oldAdaptForest;
384:   DM             oldAdapt;

390:   forest   = (DM_Forest*) dm->data;
391:   DMForestGetAdaptivityForest(dm,&oldAdapt);
392:   if (adapt != NULL && dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
393:   adaptForest    = (DM_Forest*) (adapt ? adapt->data : NULL);
394:   oldAdaptForest = (DM_Forest*) (oldAdapt ? oldAdapt->data : NULL);
395:   if (adaptForest != oldAdaptForest) {
396:     PetscSFDestroy(&forest->preCoarseToFine);
397:     PetscSFDestroy(&forest->coarseToPreFine);
398:     if (forest->clearadaptivityforest) {(forest->clearadaptivityforest)(dm);}
399:   }
400:   switch (forest->adaptPurpose) {
401:   case DM_ADAPT_DETERMINE:
402:     PetscObjectReference((PetscObject)adapt);
403:     DMDestroy(&(forest->adapt));
404:     forest->adapt = adapt;
405:     break;
406:   case DM_ADAPT_REFINE:
407:     DMSetCoarseDM(dm,adapt);
408:     break;
409:   case DM_ADAPT_COARSEN:
410:     DMSetFineDM(dm,adapt);
411:     break;
412:   default:
413:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
414:   }
415:   return(0);
416: }

418: /*@
419:   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.

421:   Not collective

423:   Input Parameter:
424: . dm - the forest

426:   Output Parameter:
427: . adapt - the forest from which dm is/was adapted

429:   Level: intermediate

431: .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose()
432: @*/
433: PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
434: {
435:   DM_Forest      *forest;

440:   forest = (DM_Forest*) dm->data;
441:   switch (forest->adaptPurpose) {
442:   case DM_ADAPT_DETERMINE:
443:     *adapt = forest->adapt;
444:     break;
445:   case DM_ADAPT_REFINE:
446:     DMGetCoarseDM(dm,adapt);
447:     break;
448:   case DM_ADAPT_COARSEN:
449:     DMGetFineDM(dm,adapt);
450:     break;
451:   default:
452:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
453:   }
454:   return(0);
455: }

457: /*@
458:   DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its
459:   source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening
460:   (DM_ADAPT_COARSEN), or undefined (DM_ADAPT_DETERMINE).  This only matters for the purposes of reference counting:
461:   during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse
462:   relationship (see DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does
463:   not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can
464:   cause a memory leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.

466:   Logically collective on dm

468:   Input Parameters:
469: + dm - the forest
470: - purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN)

472:   Level: advanced

474: .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
475: @*/
476: PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose)
477: {
478:   DM_Forest      *forest;

482:   forest = (DM_Forest*) dm->data;
483:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
484:   if (purpose != forest->adaptPurpose) {
485:     DM adapt;

487:     DMForestGetAdaptivityForest(dm,&adapt);
488:     PetscObjectReference((PetscObject)adapt);
489:     DMForestSetAdaptivityForest(dm,NULL);

491:     forest->adaptPurpose = purpose;

493:     DMForestSetAdaptivityForest(dm,adapt);
494:     DMDestroy(&adapt);
495:   }
496:   return(0);
497: }

499: /*@
500:   DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with
501:   DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening (DM_ADAPT_COARSEN), or
502:   undefined (DM_ADAPT_DETERMINE).  This only matters for the purposes of reference counting: during DMDestroy(), cyclic
503:   references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see
504:   DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does not maintain a
505:   reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory
506:   leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.

508:   Not collective

510:   Input Parameter:
511: . dm - the forest

513:   Output Parameter:
514: . purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN)

516:   Level: advanced

518: .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
519: @*/
520: PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose)
521: {
522:   DM_Forest *forest;

525:   forest   = (DM_Forest*) dm->data;
526:   *purpose = forest->adaptPurpose;
527:   return(0);
528: }

530: /*@
531:   DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
532:   cell adjacency (for the purposes of partitioning and overlap).

534:   Logically collective on dm

536:   Input Parameters:
537: + dm - the forest
538: - adjDim - default 0 (i.e., vertices determine adjacency)

540:   Level: intermediate

542: .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap()
543: @*/
544: PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
545: {
546:   PetscInt       dim;
547:   DM_Forest      *forest = (DM_Forest*) dm->data;

552:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup");
553:   if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim);
554:   DMGetDimension(dm,&dim);
555:   if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim);
556:   forest->adjDim = adjDim;
557:   return(0);
558: }

560: /*@
561:   DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that,
562:   e.g., adjacency based on facets can be specified by codimension 1 in all cases)

564:   Logically collective on dm

566:   Input Parameters:
567: + dm - the forest
568: - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices

570:   Level: intermediate

572: .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension()
573: @*/
574: PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
575: {
576:   PetscInt       dim;

581:   DMGetDimension(dm,&dim);
582:   DMForestSetAdjacencyDimension(dm,dim-adjCodim);
583:   return(0);
584: }

586: /*@
587:   DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
588:   purposes of partitioning and overlap).

590:   Not collective

592:   Input Parameter:
593: . dm - the forest

595:   Output Parameter:
596: . adjDim - default 0 (i.e., vertices determine adjacency)

598:   Level: intermediate

600: .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap()
601: @*/
602: PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
603: {
604:   DM_Forest *forest = (DM_Forest*) dm->data;

609:   *adjDim = forest->adjDim;
610:   return(0);
611: }

613: /*@
614:   DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that,
615:   e.g., adjacency based on facets can be specified by codimension 1 in all cases)

617:   Not collective

619:   Input Parameter:
620: . dm - the forest

622:   Output Parameter:
623: . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices

625:   Level: intermediate

627: .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension()
628: @*/
629: PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
630: {
631:   DM_Forest      *forest = (DM_Forest*) dm->data;
632:   PetscInt       dim;

638:   DMGetDimension(dm,&dim);
639:   *adjCodim = dim - forest->adjDim;
640:   return(0);
641: }

643: /*@
644:   DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel
645:   partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding
646:   adjacent cells

648:   Logically collective on dm

650:   Input Parameters:
651: + dm - the forest
652: - overlap - default 0

654:   Level: intermediate

656: .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
657: @*/
658: PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap)
659: {
660:   DM_Forest *forest = (DM_Forest*) dm->data;

664:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup");
665:   if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap);
666:   forest->overlap = overlap;
667:   return(0);
668: }

670: /*@
671:   DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values
672:   > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells

674:   Not collective

676:   Input Parameter:
677: . dm - the forest

679:   Output Parameter:
680: . overlap - default 0

682:   Level: intermediate

684: .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
685: @*/
686: PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap)
687: {
688:   DM_Forest *forest = (DM_Forest*) dm->data;

693:   *overlap = forest->overlap;
694:   return(0);
695: }

697: /*@
698:   DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base
699:   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest
700:   (see DMForestGetAdaptivityForest()) this limits the amount of coarsening.

702:   Logically collective on dm

704:   Input Parameters:
705: + dm - the forest
706: - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)

708:   Level: intermediate

710: .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
711: @*/
712: PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement)
713: {
714:   DM_Forest *forest = (DM_Forest*) dm->data;

718:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup");
719:   forest->minRefinement = minRefinement;
720:   return(0);
721: }

723: /*@
724:   DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see
725:   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest (see
726:   DMForestGetAdaptivityForest()), this limits the amount of coarsening.

728:   Not collective

730:   Input Parameter:
731: . dm - the forest

733:   Output Parameter:
734: . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)

736:   Level: intermediate

738: .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
739: @*/
740: PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement)
741: {
742:   DM_Forest *forest = (DM_Forest*) dm->data;

747:   *minRefinement = forest->minRefinement;
748:   return(0);
749: }

751: /*@
752:   DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base
753:   DM, see DMForestGetBaseDM()) allowed in the forest.

755:   Logically collective on dm

757:   Input Parameters:
758: + dm - the forest
759: - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)

761:   Level: intermediate

763: .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
764: @*/
765: PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement)
766: {
767:   DM_Forest *forest = (DM_Forest*) dm->data;

771:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup");
772:   forest->initRefinement = initRefinement;
773:   return(0);
774: }

776: /*@
777:   DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see
778:   DMForestGetBaseDM()) allowed in the forest.

780:   Not collective

782:   Input Parameter:
783: . dm - the forest

785:   Output Paramater:
786: . initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)

788:   Level: intermediate

790: .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
791: @*/
792: PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement)
793: {
794:   DM_Forest *forest = (DM_Forest*) dm->data;

799:   *initRefinement = forest->initRefinement;
800:   return(0);
801: }

803: /*@
804:   DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base
805:   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest
806:   (see DMForestGetAdaptivityForest()), this limits the amount of refinement.

808:   Logically collective on dm

810:   Input Parameters:
811: + dm - the forest
812: - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)

814:   Level: intermediate

816: .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM()
817: @*/
818: PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement)
819: {
820:   DM_Forest *forest = (DM_Forest*) dm->data;

824:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup");
825:   forest->maxRefinement = maxRefinement;
826:   return(0);
827: }

829: /*@
830:   DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see
831:   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest (see
832:   DMForestGetAdaptivityForest()), this limits the amount of refinement.

834:   Not collective

836:   Input Parameter:
837: . dm - the forest

839:   Output Parameter:
840: . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)

842:   Level: intermediate

844: .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
845: @*/
846: PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement)
847: {
848:   DM_Forest *forest = (DM_Forest*) dm->data;

853:   *maxRefinement = forest->maxRefinement;
854:   return(0);
855: }

857: /*@C
858:   DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes.
859:   Subtypes of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree
860:   for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to
861:   specify refinement/coarsening.

863:   Logically collective on dm

865:   Input Parameters:
866: + dm - the forest
867: - adaptStrategy - default DMFORESTADAPTALL

869:   Level: advanced

871: .seealso: DMForestGetAdaptivityStrategy()
872: @*/
873: PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy)
874: {
875:   DM_Forest      *forest = (DM_Forest*) dm->data;

880:   PetscFree(forest->adaptStrategy);
881:   PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy);
882:   return(0);
883: }

885: /*@C
886:   DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes.  Subtypes
887:   of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all
888:   processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only
889:   one process needs to specify refinement/coarsening.

891:   Not collective

893:   Input Parameter:
894: . dm - the forest

896:   Output Parameter:
897: . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL)

899:   Level: advanced

901: .seealso: DMForestSetAdaptivityStrategy()
902: @*/
903: PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
904: {
905:   DM_Forest *forest = (DM_Forest*) dm->data;

910:   *adaptStrategy = forest->adaptStrategy;
911:   return(0);
912: }

914: /*@
915:   DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning,
916:   etc.) was successful.  PETSC_FALSE indicates that the post-adaptation forest is the same as the pre-adpatation
917:   forest.  A requested adaptation may have been unsuccessful if, for example, the requested refinement would have
918:   exceeded the maximum refinement level.

920:   Collective on dm

922:   Input Parameter:

924: . dm - the post-adaptation forest

926:   Output Parameter:

928: . success - PETSC_TRUE if the post-adaptation forest is different from the pre-adaptation forest.

930:   Level: intermediate

932: .see
933: @*/
934: PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success)
935: {
936:   DM_Forest      *forest;

941:   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() has not been called yet.");
942:   forest = (DM_Forest *) dm->data;
943:   (forest->getadaptivitysuccess)(dm,success);
944:   return(0);
945: }

947: /*@
948:   DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed
949:   relating the cells of the pre-adaptation forest to the post-adaptiation forest.  After DMSetUp() is called, these transfer PetscSFs can be accessed with DMForestGetAdaptivitySF().

951:   Logically collective on dm

953:   Input Parameters:
954: + dm - the post-adaptation forest
955: - computeSF - default PETSC_TRUE

957:   Level: advanced

959: .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
960: @*/
961: PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
962: {
963:   DM_Forest *forest;

967:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called");
968:   forest                 = (DM_Forest*) dm->data;
969:   forest->computeAdaptSF = computeSF;
970:   return(0);
971: }

973: PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time)
974: {
975:   DM_Forest      *forest;

983:   forest = (DM_Forest *) dmIn->data;
984:   if (!forest->transfervec) SETERRQ(PetscObjectComm((PetscObject)dmIn),PETSC_ERR_SUP,"DMForestTransferVec() not implemented");
985:   (forest->transfervec)(dmIn,vecIn,dmOut,vecOut,useBCs,time);
986:   return(0);
987: }

989: /*@
990:   DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the
991:   pre-adaptation forest to the post-adaptiation forest.  After DMSetUp() is called, these transfer PetscSFs can be
992:   accessed with DMForestGetAdaptivitySF().

994:   Not collective

996:   Input Parameter:
997: . dm - the post-adaptation forest

999:   Output Parameter:
1000: . computeSF - default PETSC_TRUE

1002:   Level: advanced

1004: .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
1005: @*/
1006: PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1007: {
1008:   DM_Forest *forest;

1012:   forest     = (DM_Forest*) dm->data;
1013:   *computeSF = forest->computeAdaptSF;
1014:   return(0);
1015: }

1017: /*@
1018:   DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest.
1019:   Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be
1020:   some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa.  Therefore there are two
1021:   PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates
1022:   pre-adaptation fine cells to post-adaptation coarse cells.

1024:   Not collective

1026:   Input Parameter:
1027:   dm - the post-adaptation forest

1029:   Output Parameter:
1030:   preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
1031:   coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-

1033:   Level: advanced

1035: .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF()
1036: @*/
1037: PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1038: {
1039:   DM_Forest      *forest;

1044:   DMSetUp(dm);
1045:   forest = (DM_Forest*) dm->data;
1046:   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1047:   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
1048:   return(0);
1049: }

1051: /*@
1052:   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
1053:   indicate that the diameter of neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may
1054:   only support one particular choice of grading factor.

1056:   Logically collective on dm

1058:   Input Parameters:
1059: + dm - the forest
1060: - grade - the grading factor

1062:   Level: advanced

1064: .seealso: DMForestGetGradeFactor()
1065: @*/
1066: PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1067: {
1068:   DM_Forest *forest = (DM_Forest*) dm->data;

1072:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup");
1073:   forest->gradeFactor = grade;
1074:   return(0);
1075: }

1077: /*@
1078:   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
1079:   neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may only support one particular
1080:   choice of grading factor.

1082:   Not collective

1084:   Input Parameter:
1085: . dm - the forest

1087:   Output Parameter:
1088: . grade - the grading factor

1090:   Level: advanced

1092: .seealso: DMForestSetGradeFactor()
1093: @*/
1094: PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1095: {
1096:   DM_Forest *forest = (DM_Forest*) dm->data;

1101:   *grade = forest->gradeFactor;
1102:   return(0);
1103: }

1105: /*@
1106:   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
1107:   the cell weight (see DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be
1108:   (cellWeight) * (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on
1109:   its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
1110:   computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.

1112:   Logically collective on dm

1114:   Input Parameters:
1115: + dm - the forest
1116: - weightsFactors - default 1.

1118:   Level: advanced

1120: .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights()
1121: @*/
1122: PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
1123: {
1124:   DM_Forest *forest = (DM_Forest*) dm->data;

1128:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup");
1129:   forest->weightsFactor = weightsFactor;
1130:   return(0);
1131: }

1133: /*@
1134:   DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
1135:   DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be (cellWeight) *
1136:   (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on its level; a
1137:   factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation
1138:   associated with a cell is multiplied by a factor of 2 for each additional level of refinement.

1140:   Not collective

1142:   Input Parameter:
1143: . dm - the forest

1145:   Output Parameter:
1146: . weightsFactors - default 1.

1148:   Level: advanced

1150: .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights()
1151: @*/
1152: PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1153: {
1154:   DM_Forest *forest = (DM_Forest*) dm->data;

1159:   *weightsFactor = forest->weightsFactor;
1160:   return(0);
1161: }

1163: /*@
1164:   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process

1166:   Not collective

1168:   Input Parameter:
1169: . dm - the forest

1171:   Output Parameters:
1172: + cStart - the first cell on this process
1173: - cEnd - one after the final cell on this process

1175:   Level: intermediate

1177: .seealso: DMForestGetCellSF()
1178: @*/
1179: PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1180: {
1181:   DM_Forest      *forest = (DM_Forest*) dm->data;

1188:   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) {
1189:     forest->createcellchart(dm,&forest->cStart,&forest->cEnd);
1190:   }
1191:   *cStart =  forest->cStart;
1192:   *cEnd   =  forest->cEnd;
1193:   return(0);
1194: }

1196: /*@
1197:   DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes

1199:   Not collective

1201:   Input Parameter:
1202: . dm - the forest

1204:   Output Parameter:
1205: . cellSF - the PetscSF

1207:   Level: intermediate

1209: .seealso: DMForestGetCellChart()
1210: @*/
1211: PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1212: {
1213:   DM_Forest      *forest = (DM_Forest*) dm->data;

1219:   if ((!forest->cellSF) && forest->createcellsf) {
1220:     forest->createcellsf(dm,&forest->cellSF);
1221:   }
1222:   *cellSF = forest->cellSF;
1223:   return(0);
1224: }

1226: /*@C
1227:   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
1228:   DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination).  The
1229:   interpretation of the label values is up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP,
1230:   DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have been reserved as choices that should be accepted by all subtypes.

1232:   Logically collective on dm

1234:   Input Parameters:
1235: - dm - the forest
1236: + adaptLabel - the label in the pre-adaptation forest

1238:   Level: intermediate

1240: .seealso DMForestGetAdaptivityLabel()
1241: @*/
1242: PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel)
1243: {
1244:   DM_Forest      *forest = (DM_Forest*) dm->data;

1249:   adaptLabel->refct++;
1250:   if (forest->adaptLabel) {DMLabelDestroy(&forest->adaptLabel);}
1251:   forest->adaptLabel = adaptLabel;
1252:   return(0);
1253: }

1255: /*@C
1256:   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that
1257:   holds the adaptation flags (refinement, coarsening, or some combination).  The interpretation of the label values is
1258:   up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have
1259:   been reserved as choices that should be accepted by all subtypes.

1261:   Not collective

1263:   Input Parameter:
1264: . dm - the forest

1266:   Output Parameter:
1267: . adaptLabel - the name of the label in the pre-adaptation forest

1269:   Level: intermediate

1271: .seealso DMForestSetAdaptivityLabel()
1272: @*/
1273: PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel)
1274: {
1275:   DM_Forest *forest = (DM_Forest*) dm->data;

1279:   *adaptLabel = forest->adaptLabel;
1280:   return(0);
1281: }

1283: /*@
1284:   DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
1285:   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
1286:   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
1287:   of 1.

1289:   Logically collective on dm

1291:   Input Parameters:
1292: + dm - the forest
1293: . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
1294: - copyMode - how weights should reference weights

1296:   Level: advanced

1298: .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity()
1299: @*/
1300: PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1301: {
1302:   DM_Forest      *forest = (DM_Forest*) dm->data;
1303:   PetscInt       cStart, cEnd;

1308:   DMForestGetCellChart(dm,&cStart,&cEnd);
1309:   if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd);
1310:   if (copyMode == PETSC_COPY_VALUES) {
1311:     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) {
1312:       PetscMalloc1(cEnd-cStart,&forest->cellWeights);
1313:     }
1314:     PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));
1315:     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
1316:     return(0);
1317:   }
1318:   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) {
1319:     PetscFree(forest->cellWeights);
1320:   }
1321:   forest->cellWeights         = weights;
1322:   forest->cellWeightsCopyMode = copyMode;
1323:   return(0);
1324: }

1326: /*@
1327:   DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
1328:   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
1329:   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
1330:   of 1.

1332:   Not collective

1334:   Input Parameter:
1335: . dm - the forest

1337:   Output Parameter:
1338: . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.

1340:   Level: advanced

1342: .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity()
1343: @*/
1344: PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1345: {
1346:   DM_Forest *forest = (DM_Forest*) dm->data;

1351:   *weights = forest->cellWeights;
1352:   return(0);
1353: }

1355: /*@
1356:   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
1357:   a pre-adaptation forest (see DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each
1358:   process's cells to the process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that
1359:   the current process should not have any cells after repartitioning.

1361:   Logically Collective on dm

1363:   Input parameters:
1364: + dm - the forest
1365: - capacity - this process's capacity

1367:   Level: advanced

1369: .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
1370: @*/
1371: PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1372: {
1373:   DM_Forest *forest = (DM_Forest*) dm->data;

1377:   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup");
1378:   if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity);
1379:   forest->weightCapacity = capacity;
1380:   return(0);
1381: }

1383: /*@
1384:   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
1385:   DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each process's cells to the
1386:   process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that the current process
1387:   should not have any cells after repartitioning.

1389:   Not collective

1391:   Input parameter:
1392: . dm - the forest

1394:   Output parameter:
1395: . capacity - this process's capacity

1397:   Level: advanced

1399: .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
1400: @*/
1401: PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1402: {
1403:   DM_Forest *forest = (DM_Forest*) dm->data;

1408:   *capacity = forest->weightCapacity;
1409:   return(0);
1410: }

1412: PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm)
1413: {
1414:   DM_Forest                  *forest = (DM_Forest*) dm->data;
1415:   PetscBool                  flg, flg1, flg2, flg3, flg4;
1416:   DMForestTopology           oldTopo;
1417:   char                       stringBuffer[256];
1418:   PetscViewer                viewer;
1419:   PetscViewerFormat          format;
1420:   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1421:   PetscReal                  weightsFactor;
1422:   DMForestAdaptivityStrategy adaptStrategy;
1423:   PetscErrorCode             ierr;

1427:   forest->setfromoptionscalled = PETSC_TRUE;
1428:   DMForestGetTopology(dm, &oldTopo);
1429:   PetscOptionsHead(PetscOptionsObject,"DMForest Options");
1430:   PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);
1431:   PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);
1432:   PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);
1433:   PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);
1434:   if ((PetscInt) flg1 + (PetscInt) flg2 + (PetscInt) flg3 + (PetscInt) flg4 > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Specify only one of -dm_forest_{topology,base_dm,coarse_forest,fine_forest}");
1435:   if (flg1) {
1436:     DMForestSetTopology(dm,(DMForestTopology)stringBuffer);
1437:     DMForestSetBaseDM(dm,NULL);
1438:     DMForestSetAdaptivityForest(dm,NULL);
1439:   }
1440:   if (flg2) {
1441:     DM base;

1443:     DMCreate(PetscObjectComm((PetscObject)dm),&base);
1444:     PetscViewerPushFormat(viewer,format);
1445:     DMLoad(base,viewer);
1446:     PetscViewerDestroy(&viewer);
1447:     DMForestSetBaseDM(dm,base);
1448:     DMDestroy(&base);
1449:     DMForestSetTopology(dm,NULL);
1450:     DMForestSetAdaptivityForest(dm,NULL);
1451:   }
1452:   if (flg3) {
1453:     DM coarse;

1455:     DMCreate(PetscObjectComm((PetscObject)dm),&coarse);
1456:     PetscViewerPushFormat(viewer,format);
1457:     DMLoad(coarse,viewer);
1458:     PetscViewerDestroy(&viewer);
1459:     DMForestSetAdaptivityForest(dm,coarse);
1460:     DMDestroy(&coarse);
1461:     DMForestSetTopology(dm,NULL);
1462:     DMForestSetBaseDM(dm,NULL);
1463:   }
1464:   if (flg4) {
1465:     DM fine;

1467:     DMCreate(PetscObjectComm((PetscObject)dm),&fine);
1468:     PetscViewerPushFormat(viewer,format);
1469:     DMLoad(fine,viewer);
1470:     PetscViewerDestroy(&viewer);
1471:     DMForestSetAdaptivityForest(dm,fine);
1472:     DMDestroy(&fine);
1473:     DMForestSetTopology(dm,NULL);
1474:     DMForestSetBaseDM(dm,NULL);
1475:   }
1476:   DMForestGetAdjacencyDimension(dm,&adjDim);
1477:   PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);
1478:   if (flg) {
1479:     DMForestSetAdjacencyDimension(dm,adjDim);
1480:   } else {
1481:     DMForestGetAdjacencyCodimension(dm,&adjCodim);
1482:     PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);
1483:     if (flg) {
1484:       DMForestSetAdjacencyCodimension(dm,adjCodim);
1485:     }
1486:   }
1487:   DMForestGetPartitionOverlap(dm,&overlap);
1488:   PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);
1489:   if (flg) {
1490:     DMForestSetPartitionOverlap(dm,overlap);
1491:   }
1492: #if 0
1493:   PetscOptionsInt("-dm_refine","equivalent to -dm_forest_set_minimum_refinement and -dm_forest_set_initial_refinement with the same value",NULL,minRefinement,&minRefinement,&flg);
1494:   if (flg) {
1495:     DMForestSetMinimumRefinement(dm,minRefinement);
1496:     DMForestSetInitialRefinement(dm,minRefinement);
1497:   }
1498:   PetscOptionsInt("-dm_refine_hierarchy","equivalent to -dm_forest_set_minimum_refinement 0 and -dm_forest_set_initial_refinement",NULL,initRefinement,&initRefinement,&flg);
1499:   if (flg) {
1500:     DMForestSetMinimumRefinement(dm,0);
1501:     DMForestSetInitialRefinement(dm,initRefinement);
1502:   }
1503: #endif
1504:   DMForestGetMinimumRefinement(dm,&minRefinement);
1505:   PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);
1506:   if (flg) {
1507:     DMForestSetMinimumRefinement(dm,minRefinement);
1508:   }
1509:   DMForestGetInitialRefinement(dm,&initRefinement);
1510:   PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);
1511:   if (flg) {
1512:     DMForestSetInitialRefinement(dm,initRefinement);
1513:   }
1514:   DMForestGetMaximumRefinement(dm,&maxRefinement);
1515:   PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);
1516:   if (flg) {
1517:     DMForestSetMaximumRefinement(dm,maxRefinement);
1518:   }
1519:   DMForestGetAdaptivityStrategy(dm,&adaptStrategy);
1520:   PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);
1521:   if (flg) {
1522:     DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);
1523:   }
1524:   DMForestGetGradeFactor(dm,&grade);
1525:   PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);
1526:   if (flg) {
1527:     DMForestSetGradeFactor(dm,grade);
1528:   }
1529:   DMForestGetCellWeightFactor(dm,&weightsFactor);
1530:   PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);
1531:   if (flg) {
1532:     DMForestSetCellWeightFactor(dm,weightsFactor);
1533:   }
1534:   PetscOptionsTail();
1535:   return(0);
1536: }

1538: PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
1539: {

1543:   if (subdm) {DMClone(dm, subdm);}
1544:   DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);
1545:   return(0);
1546: }

1548: PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
1549: {
1550:   DMLabel        refine;
1551:   DM             fineDM;

1555:   DMGetFineDM(dm,&fineDM);
1556:   if (fineDM) {
1557:     PetscObjectReference((PetscObject)fineDM);
1558:     *dmRefined = fineDM;
1559:     return(0);
1560:   }
1561:   DMForestTemplate(dm,comm,dmRefined);
1562:   DMGetLabel(dm,"refine",&refine);
1563:   if (!refine) {
1564:     DMLabelCreate("refine",&refine);
1565:     DMLabelSetDefaultValue(refine,DM_ADAPT_REFINE);
1566:   }
1567:   else {
1568:     refine->refct++;
1569:   }
1570:   DMForestSetAdaptivityLabel(*dmRefined,refine);
1571:   DMLabelDestroy(&refine);
1572:   return(0);
1573: }

1575: PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
1576: {
1577:   DMLabel        coarsen;
1578:   DM             coarseDM;

1582:   {
1583:     PetscMPIInt mpiComparison;
1584:     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);

1586:     MPI_Comm_compare(comm, dmcomm, &mpiComparison);
1587:     if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet");
1588:   }
1589:   DMGetCoarseDM(dm,&coarseDM);
1590:   if (coarseDM) {
1591:     PetscObjectReference((PetscObject)coarseDM);
1592:     *dmCoarsened = coarseDM;
1593:     return(0);
1594:   }
1595:   DMForestTemplate(dm,comm,dmCoarsened);
1596:   DMForestSetAdaptivityPurpose(coarseDM,DM_ADAPT_COARSEN);
1597:   DMGetLabel(dm,"coarsen",&coarsen);
1598:   if (!coarsen) {
1599:     DMLabelCreate("coarsen",&coarsen);
1600:     DMLabelSetDefaultValue(coarsen,DM_ADAPT_COARSEN);
1601:   } else {
1602:     coarsen->refct++;
1603:   }
1604:   DMForestSetAdaptivityLabel(*dmCoarsened,coarsen);
1605:   DMLabelDestroy(&coarsen);
1606:   return(0);
1607: }

1609: static PetscErrorCode DMAdaptLabel_Forest(DM dm, DMLabel label, DM *adaptedDM)
1610: {
1611:   PetscBool      success;

1615:   DMForestTemplate(dm,PetscObjectComm((PetscObject)dm),adaptedDM);
1616:   DMForestSetAdaptivityLabel(*adaptedDM,label);
1617:   DMSetUp(*adaptedDM);
1618:   DMForestGetAdaptivitySuccess(*adaptedDM,&success);
1619:   if (!success) {
1620:     DMDestroy(adaptedDM);
1621:     *adaptedDM = NULL;
1622:   }
1623:   return(0);
1624: }

1626: static PetscErrorCode DMInitialize_Forest(DM dm)
1627: {

1631:   PetscMemzero(dm->ops,sizeof(*(dm->ops)));

1633:   dm->ops->clone          = DMClone_Forest;
1634:   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1635:   dm->ops->destroy        = DMDestroy_Forest;
1636:   dm->ops->createsubdm    = DMCreateSubDM_Forest;
1637:   dm->ops->refine         = DMRefine_Forest;
1638:   dm->ops->coarsen        = DMCoarsen_Forest;
1639:   dm->ops->adaptlabel     = DMAdaptLabel_Forest;
1640:   return(0);
1641: }

1643: /*MC

1645:      DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh.  Forests usually have a base DM
1646:   (see DMForestGetBaseDM()), from which it is refined.  The refinement and partitioning of forests is considered
1647:   immutable after DMSetUp() is called.  To adapt a mesh, one should call DMForestTemplate() to create a new mesh that
1648:   will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new
1649:   mesh.

1651:   To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the
1652:   previous mesh whose values indicate which cells should be refined (DM_ADAPT_REFINE) or coarsened (DM_ADAPT_COARSEN)
1653:   and how (subtypes are free to allow additional values for things like anisotropic refinement).  The label should be
1654:   given to the *new* mesh with DMForestSetAdaptivityLabel().

1656:   Level: advanced

1658: .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel()
1659: M*/

1661: PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1662: {
1663:   DM_Forest      *forest;

1668:   PetscNewLog(dm,&forest);
1669:   dm->dim                      = 0;
1670:   dm->data                     = forest;
1671:   forest->refct                = 1;
1672:   forest->data                 = NULL;
1673:   forest->setfromoptionscalled = PETSC_FALSE;
1674:   forest->topology             = NULL;
1675:   forest->adapt                = NULL;
1676:   forest->base                 = NULL;
1677:   forest->adaptPurpose         = DM_ADAPT_DETERMINE;
1678:   forest->adjDim               = PETSC_DEFAULT;
1679:   forest->overlap              = PETSC_DEFAULT;
1680:   forest->minRefinement        = PETSC_DEFAULT;
1681:   forest->maxRefinement        = PETSC_DEFAULT;
1682:   forest->initRefinement       = PETSC_DEFAULT;
1683:   forest->cStart               = PETSC_DETERMINE;
1684:   forest->cEnd                 = PETSC_DETERMINE;
1685:   forest->cellSF               = NULL;
1686:   forest->adaptLabel           = NULL;
1687:   forest->gradeFactor          = 2;
1688:   forest->cellWeights          = NULL;
1689:   forest->cellWeightsCopyMode  = PETSC_USE_POINTER;
1690:   forest->weightsFactor        = 1.;
1691:   forest->weightCapacity       = 1.;
1692:   DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);
1693:   DMInitialize_Forest(dm);
1694:   return(0);
1695: }