Actual source code: forest.c

petsc-3.9.4 2018-09-11
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;
385:   PetscBool      isForest;

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

421: /*@
422:   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.

424:   Not collective

426:   Input Parameter:
427: . dm - the forest

429:   Output Parameter:
430: . adapt - the forest from which dm is/was adapted

432:   Level: intermediate

434: .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose()
435: @*/
436: PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
437: {
438:   DM_Forest      *forest;

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

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

469:   Logically collective on dm

471:   Input Parameters:
472: + dm - the forest
473: - purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN)

475:   Level: advanced

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

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

490:     DMForestGetAdaptivityForest(dm,&adapt);
491:     PetscObjectReference((PetscObject)adapt);
492:     DMForestSetAdaptivityForest(dm,NULL);

494:     forest->adaptPurpose = purpose;

496:     DMForestSetAdaptivityForest(dm,adapt);
497:     DMDestroy(&adapt);
498:   }
499:   return(0);
500: }

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

511:   Not collective

513:   Input Parameter:
514: . dm - the forest

516:   Output Parameter:
517: . purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN)

519:   Level: advanced

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

528:   forest   = (DM_Forest*) dm->data;
529:   *purpose = forest->adaptPurpose;
530:   return(0);
531: }

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

537:   Logically collective on dm

539:   Input Parameters:
540: + dm - the forest
541: - adjDim - default 0 (i.e., vertices determine adjacency)

543:   Level: intermediate

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

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

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

567:   Logically collective on dm

569:   Input Parameters:
570: + dm - the forest
571: - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices

573:   Level: intermediate

575: .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension()
576: @*/
577: PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
578: {
579:   PetscInt       dim;

584:   DMGetDimension(dm,&dim);
585:   DMForestSetAdjacencyDimension(dm,dim-adjCodim);
586:   return(0);
587: }

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

593:   Not collective

595:   Input Parameter:
596: . dm - the forest

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

601:   Level: intermediate

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

612:   *adjDim = forest->adjDim;
613:   return(0);
614: }

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

620:   Not collective

622:   Input Parameter:
623: . dm - the forest

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

628:   Level: intermediate

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

641:   DMGetDimension(dm,&dim);
642:   *adjCodim = dim - forest->adjDim;
643:   return(0);
644: }

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

651:   Logically collective on dm

653:   Input Parameters:
654: + dm - the forest
655: - overlap - default 0

657:   Level: intermediate

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

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

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

677:   Not collective

679:   Input Parameter:
680: . dm - the forest

682:   Output Parameter:
683: . overlap - default 0

685:   Level: intermediate

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

696:   *overlap = forest->overlap;
697:   return(0);
698: }

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

705:   Logically collective on dm

707:   Input Parameters:
708: + dm - the forest
709: - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)

711:   Level: intermediate

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

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

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

731:   Not collective

733:   Input Parameter:
734: . dm - the forest

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

739:   Level: intermediate

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

750:   *minRefinement = forest->minRefinement;
751:   return(0);
752: }

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

758:   Logically collective on dm

760:   Input Parameters:
761: + dm - the forest
762: - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)

764:   Level: intermediate

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

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

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

783:   Not collective

785:   Input Parameter:
786: . dm - the forest

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

791:   Level: intermediate

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

802:   *initRefinement = forest->initRefinement;
803:   return(0);
804: }

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

811:   Logically collective on dm

813:   Input Parameters:
814: + dm - the forest
815: - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)

817:   Level: intermediate

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

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

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

837:   Not collective

839:   Input Parameter:
840: . dm - the forest

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

845:   Level: intermediate

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

856:   *maxRefinement = forest->maxRefinement;
857:   return(0);
858: }

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

866:   Logically collective on dm

868:   Input Parameters:
869: + dm - the forest
870: - adaptStrategy - default DMFORESTADAPTALL

872:   Level: advanced

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

883:   PetscFree(forest->adaptStrategy);
884:   PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy);
885:   return(0);
886: }

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

894:   Not collective

896:   Input Parameter:
897: . dm - the forest

899:   Output Parameter:
900: . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL)

902:   Level: advanced

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

913:   *adaptStrategy = forest->adaptStrategy;
914:   return(0);
915: }

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

923:   Collective on dm

925:   Input Parameter:

927: . dm - the post-adaptation forest

929:   Output Parameter:

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

933:   Level: intermediate

935: .see
936: @*/
937: PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success)
938: {
939:   DM_Forest      *forest;

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

950: /*@
951:   DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed
952:   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().

954:   Logically collective on dm

956:   Input Parameters:
957: + dm - the post-adaptation forest
958: - computeSF - default PETSC_TRUE

960:   Level: advanced

962: .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
963: @*/
964: PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
965: {
966:   DM_Forest *forest;

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

976: PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time)
977: {
978:   DM_Forest      *forest;

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

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

997:   Not collective

999:   Input Parameter:
1000: . dm - the post-adaptation forest

1002:   Output Parameter:
1003: . computeSF - default PETSC_TRUE

1005:   Level: advanced

1007: .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
1008: @*/
1009: PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1010: {
1011:   DM_Forest *forest;

1015:   forest     = (DM_Forest*) dm->data;
1016:   *computeSF = forest->computeAdaptSF;
1017:   return(0);
1018: }

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

1027:   Not collective

1029:   Input Parameter:
1030:   dm - the post-adaptation forest

1032:   Output Parameter:
1033:   preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
1034:   coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-

1036:   Level: advanced

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

1047:   DMSetUp(dm);
1048:   forest = (DM_Forest*) dm->data;
1049:   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1050:   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
1051:   return(0);
1052: }

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

1059:   Logically collective on dm

1061:   Input Parameters:
1062: + dm - the forest
1063: - grade - the grading factor

1065:   Level: advanced

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

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

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

1085:   Not collective

1087:   Input Parameter:
1088: . dm - the forest

1090:   Output Parameter:
1091: . grade - the grading factor

1093:   Level: advanced

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

1104:   *grade = forest->gradeFactor;
1105:   return(0);
1106: }

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

1115:   Logically collective on dm

1117:   Input Parameters:
1118: + dm - the forest
1119: - weightsFactors - default 1.

1121:   Level: advanced

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

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

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

1143:   Not collective

1145:   Input Parameter:
1146: . dm - the forest

1148:   Output Parameter:
1149: . weightsFactors - default 1.

1151:   Level: advanced

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

1162:   *weightsFactor = forest->weightsFactor;
1163:   return(0);
1164: }

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

1169:   Not collective

1171:   Input Parameter:
1172: . dm - the forest

1174:   Output Parameters:
1175: + cStart - the first cell on this process
1176: - cEnd - one after the final cell on this process

1178:   Level: intermediate

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

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

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

1202:   Not collective

1204:   Input Parameter:
1205: . dm - the forest

1207:   Output Parameter:
1208: . cellSF - the PetscSF

1210:   Level: intermediate

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

1222:   if ((!forest->cellSF) && forest->createcellsf) {
1223:     forest->createcellsf(dm,&forest->cellSF);
1224:   }
1225:   *cellSF = forest->cellSF;
1226:   return(0);
1227: }

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

1235:   Logically collective on dm

1237:   Input Parameters:
1238: - dm - the forest
1239: + adaptLabel - the label in the pre-adaptation forest

1241:   Level: intermediate

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

1252:   adaptLabel->refct++;
1253:   if (forest->adaptLabel) {DMLabelDestroy(&forest->adaptLabel);}
1254:   forest->adaptLabel = adaptLabel;
1255:   return(0);
1256: }

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

1264:   Not collective

1266:   Input Parameter:
1267: . dm - the forest

1269:   Output Parameter:
1270: . adaptLabel - the name of the label in the pre-adaptation forest

1272:   Level: intermediate

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

1282:   *adaptLabel = forest->adaptLabel;
1283:   return(0);
1284: }

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

1292:   Logically collective on dm

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

1299:   Level: advanced

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

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

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

1335:   Not collective

1337:   Input Parameter:
1338: . dm - the forest

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

1343:   Level: advanced

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

1354:   *weights = forest->cellWeights;
1355:   return(0);
1356: }

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

1364:   Logically Collective on dm

1366:   Input parameters:
1367: + dm - the forest
1368: - capacity - this process's capacity

1370:   Level: advanced

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

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

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

1392:   Not collective

1394:   Input parameter:
1395: . dm - the forest

1397:   Output parameter:
1398: . capacity - this process's capacity

1400:   Level: advanced

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

1411:   *capacity = forest->weightCapacity;
1412:   return(0);
1413: }

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

1430:   forest->setfromoptionscalled = PETSC_TRUE;
1431:   DMForestGetTopology(dm, &oldTopo);
1432:   PetscOptionsHead(PetscOptionsObject,"DMForest Options");
1433:   PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);
1434:   PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);
1435:   PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);
1436:   PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);
1437:   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}");
1438:   if (flg1) {
1439:     DMForestSetTopology(dm,(DMForestTopology)stringBuffer);
1440:     DMForestSetBaseDM(dm,NULL);
1441:     DMForestSetAdaptivityForest(dm,NULL);
1442:   }
1443:   if (flg2) {
1444:     DM base;

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

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

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

1541: PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1542: {

1546:   if (subdm) {DMClone(dm, subdm);}
1547:   DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);
1548:   return(0);
1549: }

1551: PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
1552: {
1553:   DMLabel        refine;
1554:   DM             fineDM;

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

1578: PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
1579: {
1580:   DMLabel        coarsen;
1581:   DM             coarseDM;

1585:   {
1586:     PetscMPIInt mpiComparison;
1587:     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);

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

1612: static PetscErrorCode DMAdaptLabel_Forest(DM dm, DMLabel label, DM *adaptedDM)
1613: {
1614:   PetscBool      success;

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

1629: static PetscErrorCode DMInitialize_Forest(DM dm)
1630: {

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

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

1646: /*MC

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

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

1659:   Level: advanced

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

1664: PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1665: {
1666:   DM_Forest      *forest;

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