Actual source code: forest.c

petsc-3.10.5 2019-03-28
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:
285:     Currently the base DM must be a DMPLEX

287:   Level: intermediate

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

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

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

318:   Not collective

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

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

326:   Level: intermediate

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

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

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

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

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

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

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

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

372:   Logically collective on dm

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

378:   Level: intermediate

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

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

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

425:   Not collective

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

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

433:   Level: intermediate

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

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

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

470:   Logically collective on dm

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

476:   Level: advanced

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

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

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

495:     forest->adaptPurpose = purpose;

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

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

512:   Not collective

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

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

520:   Level: advanced

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

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

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

538:   Logically collective on dm

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

544:   Level: intermediate

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

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

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

568:   Logically collective on dm

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

574:   Level: intermediate

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

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

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

594:   Not collective

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

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

602:   Level: intermediate

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

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

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

621:   Not collective

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

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

629:   Level: intermediate

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

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

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

652:   Logically collective on dm

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

658:   Level: intermediate

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

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

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

678:   Not collective

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

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

686:   Level: intermediate

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

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

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

706:   Logically collective on dm

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

712:   Level: intermediate

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

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

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

732:   Not collective

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

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

740:   Level: intermediate

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

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

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

759:   Logically collective on dm

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

765:   Level: intermediate

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

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

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

784:   Not collective

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

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

792:   Level: intermediate

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

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

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

812:   Logically collective on dm

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

818:   Level: intermediate

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

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

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

838:   Not collective

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

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

846:   Level: intermediate

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

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

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

867:   Logically collective on dm

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

873:   Level: advanced

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

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

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

895:   Not collective

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

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

903:   Level: advanced

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

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

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

924:   Collective on dm

926:   Input Parameter:

928: . dm - the post-adaptation forest

930:   Output Parameter:

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

934:   Level: intermediate

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

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

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

955:   Logically collective on dm

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

961:   Level: advanced

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

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

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

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

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

998:   Not collective

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

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

1006:   Level: advanced

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

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

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

1028:   Not collective

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

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

1037:   Level: advanced

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

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

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

1060:   Logically collective on dm

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

1066:   Level: advanced

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

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

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

1086:   Not collective

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

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

1094:   Level: advanced

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

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

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

1116:   Logically collective on dm

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

1122:   Level: advanced

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

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

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

1144:   Not collective

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

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

1152:   Level: advanced

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

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

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

1170:   Not collective

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

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

1179:   Level: intermediate

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

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

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

1203:   Not collective

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

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

1211:   Level: intermediate

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

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

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

1236:   Logically collective on dm

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

1242:   Level: intermediate

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

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

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

1265:   Not collective

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

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

1273:   Level: intermediate

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

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

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

1293:   Logically collective on dm

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

1300:   Level: advanced

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

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

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

1336:   Not collective

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

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

1344:   Level: advanced

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

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

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

1365:   Logically Collective on dm

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

1371:   Level: advanced

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

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

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

1393:   Not collective

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

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

1401:   Level: advanced

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1647: /*MC

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

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

1660:   Level: advanced

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

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

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