xref: /petsc/src/dm/impls/forest/forest.c (revision dce8aeba1c9b69b19f651c53d8a6b674bd7e9cbd)
19be51f97SToby Isaac #include <petsc/private/dmforestimpl.h> /*I "petscdmforest.h" I*/
29be51f97SToby Isaac #include <petsc/private/dmimpl.h>       /*I "petscdm.h" I*/
3a1b0c543SToby Isaac #include <petsc/private/dmlabelimpl.h>  /*I "petscdmlabel.h" I*/
4ef19d27cSToby Isaac #include <petscsf.h>
5db4d5e8cSToby Isaac 
627d4645fSToby Isaac PetscBool DMForestPackageInitialized = PETSC_FALSE;
727d4645fSToby Isaac 
827d4645fSToby Isaac typedef struct _DMForestTypeLink *DMForestTypeLink;
927d4645fSToby Isaac 
109371c9d4SSatish Balay struct _DMForestTypeLink {
1127d4645fSToby Isaac   char            *name;
1227d4645fSToby Isaac   DMForestTypeLink next;
1327d4645fSToby Isaac };
1427d4645fSToby Isaac 
1527d4645fSToby Isaac DMForestTypeLink DMForestTypeList;
1627d4645fSToby Isaac 
17d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMForestPackageFinalize(void)
18d71ae5a4SJacob Faibussowitsch {
1927d4645fSToby Isaac   DMForestTypeLink oldLink, link = DMForestTypeList;
2027d4645fSToby Isaac 
2127d4645fSToby Isaac   PetscFunctionBegin;
2227d4645fSToby Isaac   while (link) {
2327d4645fSToby Isaac     oldLink = link;
249566063dSJacob Faibussowitsch     PetscCall(PetscFree(oldLink->name));
2527d4645fSToby Isaac     link = oldLink->next;
269566063dSJacob Faibussowitsch     PetscCall(PetscFree(oldLink));
2727d4645fSToby Isaac   }
2827d4645fSToby Isaac   PetscFunctionReturn(0);
2927d4645fSToby Isaac }
3027d4645fSToby Isaac 
31d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMForestPackageInitialize(void)
32d71ae5a4SJacob Faibussowitsch {
3327d4645fSToby Isaac   PetscFunctionBegin;
3427d4645fSToby Isaac   if (DMForestPackageInitialized) PetscFunctionReturn(0);
3527d4645fSToby Isaac   DMForestPackageInitialized = PETSC_TRUE;
36f885a11aSToby Isaac 
379566063dSJacob Faibussowitsch   PetscCall(DMForestRegisterType(DMFOREST));
389566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(DMForestPackageFinalize));
3927d4645fSToby Isaac   PetscFunctionReturn(0);
4027d4645fSToby Isaac }
4127d4645fSToby Isaac 
429be51f97SToby Isaac /*@C
43*dce8aebaSBarry Smith   DMForestRegisterType - Registers a `DMType` as a subtype of `DMFOREST` (so that `DMIsForest()` will be correct)
449be51f97SToby Isaac 
459be51f97SToby Isaac   Not Collective
469be51f97SToby Isaac 
479be51f97SToby Isaac   Input parameter:
489be51f97SToby Isaac . name - the name of the type
499be51f97SToby Isaac 
509be51f97SToby Isaac   Level: advanced
519be51f97SToby Isaac 
52db781477SPatrick Sanan .seealso: `DMFOREST`, `DMIsForest()`
539be51f97SToby Isaac @*/
54d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestRegisterType(DMType name)
55d71ae5a4SJacob Faibussowitsch {
5627d4645fSToby Isaac   DMForestTypeLink link;
5727d4645fSToby Isaac 
5827d4645fSToby Isaac   PetscFunctionBegin;
599566063dSJacob Faibussowitsch   PetscCall(DMForestPackageInitialize());
609566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
619566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &link->name));
6227d4645fSToby Isaac   link->next       = DMForestTypeList;
6327d4645fSToby Isaac   DMForestTypeList = link;
6427d4645fSToby Isaac   PetscFunctionReturn(0);
6527d4645fSToby Isaac }
6627d4645fSToby Isaac 
679be51f97SToby Isaac /*@
689be51f97SToby Isaac   DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes
699be51f97SToby Isaac 
709be51f97SToby Isaac   Not Collective
719be51f97SToby Isaac 
729be51f97SToby Isaac   Input parameter:
739be51f97SToby Isaac . dm - the DM object
749be51f97SToby Isaac 
759be51f97SToby Isaac   Output parameter:
769be51f97SToby Isaac . isForest - whether dm is a subtype of DMFOREST
779be51f97SToby Isaac 
789be51f97SToby Isaac   Level: intermediate
799be51f97SToby Isaac 
80db781477SPatrick Sanan .seealso: `DMFOREST`, `DMForestRegisterType()`
819be51f97SToby Isaac @*/
82d71ae5a4SJacob Faibussowitsch PetscErrorCode DMIsForest(DM dm, PetscBool *isForest)
83d71ae5a4SJacob Faibussowitsch {
8427d4645fSToby Isaac   DMForestTypeLink link = DMForestTypeList;
8527d4645fSToby Isaac 
8627d4645fSToby Isaac   PetscFunctionBegin;
8727d4645fSToby Isaac   while (link) {
8827d4645fSToby Isaac     PetscBool sameType;
899566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)dm, link->name, &sameType));
9027d4645fSToby Isaac     if (sameType) {
9127d4645fSToby Isaac       *isForest = PETSC_TRUE;
9227d4645fSToby Isaac       PetscFunctionReturn(0);
9327d4645fSToby Isaac     }
9427d4645fSToby Isaac     link = link->next;
9527d4645fSToby Isaac   }
9627d4645fSToby Isaac   *isForest = PETSC_FALSE;
9727d4645fSToby Isaac   PetscFunctionReturn(0);
9827d4645fSToby Isaac }
9927d4645fSToby Isaac 
1009be51f97SToby Isaac /*@
1019be51f97SToby Isaac   DMForestTemplate - Create a new DM that will be adapted from a source DM.  The new DM reproduces the configuration
1029be51f97SToby Isaac   of the source, but is not yet setup, so that the user can then define only the ways that the new DM should differ
1039be51f97SToby Isaac   (by, e.g., refinement or repartitioning).  The source DM is also set as the adaptivity source DM of the new DM (see
1049be51f97SToby Isaac   DMForestSetAdaptivityForest()).
1059be51f97SToby Isaac 
1069be51f97SToby Isaac   Collective on dm
1079be51f97SToby Isaac 
1089be51f97SToby Isaac   Input Parameters:
1099be51f97SToby Isaac + dm - the source DM object
1109be51f97SToby Isaac - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen())
1119be51f97SToby Isaac 
1129be51f97SToby Isaac   Output Parameter:
1139be51f97SToby Isaac . tdm - the new DM object
1149be51f97SToby Isaac 
1159be51f97SToby Isaac   Level: intermediate
1169be51f97SToby Isaac 
117db781477SPatrick Sanan .seealso: `DMForestSetAdaptivityForest()`
1189be51f97SToby Isaac @*/
119d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm)
120d71ae5a4SJacob Faibussowitsch {
121a0452a8eSToby Isaac   DM_Forest                 *forest = (DM_Forest *)dm->data;
12220e8089bSToby Isaac   DMType                     type;
123a0452a8eSToby Isaac   DM                         base;
124a0452a8eSToby Isaac   DMForestTopology           topology;
12505e99e11SStefano Zampini   MatType                    mtype;
126a0452a8eSToby Isaac   PetscInt                   dim, overlap, ref, factor;
127a0452a8eSToby Isaac   DMForestAdaptivityStrategy strat;
128795844e7SToby Isaac   void                      *ctx;
12949fc9a2fSToby Isaac   PetscErrorCode (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *);
1303e58adeeSToby Isaac   void *mapCtx;
131a0452a8eSToby Isaac 
132a0452a8eSToby Isaac   PetscFunctionBegin;
133a0452a8eSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1349566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), tdm));
1359566063dSJacob Faibussowitsch   PetscCall(DMGetType(dm, &type));
1369566063dSJacob Faibussowitsch   PetscCall(DMSetType(*tdm, type));
1379566063dSJacob Faibussowitsch   PetscCall(DMForestGetBaseDM(dm, &base));
1389566063dSJacob Faibussowitsch   PetscCall(DMForestSetBaseDM(*tdm, base));
1399566063dSJacob Faibussowitsch   PetscCall(DMForestGetTopology(dm, &topology));
1409566063dSJacob Faibussowitsch   PetscCall(DMForestSetTopology(*tdm, topology));
1419566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdjacencyDimension(dm, &dim));
1429566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdjacencyDimension(*tdm, dim));
1439566063dSJacob Faibussowitsch   PetscCall(DMForestGetPartitionOverlap(dm, &overlap));
1449566063dSJacob Faibussowitsch   PetscCall(DMForestSetPartitionOverlap(*tdm, overlap));
1459566063dSJacob Faibussowitsch   PetscCall(DMForestGetMinimumRefinement(dm, &ref));
1469566063dSJacob Faibussowitsch   PetscCall(DMForestSetMinimumRefinement(*tdm, ref));
1479566063dSJacob Faibussowitsch   PetscCall(DMForestGetMaximumRefinement(dm, &ref));
1489566063dSJacob Faibussowitsch   PetscCall(DMForestSetMaximumRefinement(*tdm, ref));
1499566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivityStrategy(dm, &strat));
1509566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityStrategy(*tdm, strat));
1519566063dSJacob Faibussowitsch   PetscCall(DMForestGetGradeFactor(dm, &factor));
1529566063dSJacob Faibussowitsch   PetscCall(DMForestSetGradeFactor(*tdm, factor));
1539566063dSJacob Faibussowitsch   PetscCall(DMForestGetBaseCoordinateMapping(dm, &map, &mapCtx));
1549566063dSJacob Faibussowitsch   PetscCall(DMForestSetBaseCoordinateMapping(*tdm, map, mapCtx));
1551baa6e33SBarry Smith   if (forest->ftemplate) PetscCall((*forest->ftemplate)(dm, *tdm));
1569566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityForest(*tdm, dm));
1579566063dSJacob Faibussowitsch   PetscCall(DMCopyDisc(dm, *tdm));
1589566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dm, &ctx));
1599566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*tdm, &ctx));
16090b157c4SStefano Zampini   {
1614fb89dddSMatthew G. Knepley     const PetscReal *maxCell, *L, *Lstart;
162795844e7SToby Isaac 
1634fb89dddSMatthew G. Knepley     PetscCall(DMGetPeriodicity(dm, &maxCell, &Lstart, &L));
1644fb89dddSMatthew G. Knepley     PetscCall(DMSetPeriodicity(*tdm, maxCell, Lstart, L));
165795844e7SToby Isaac   }
1669566063dSJacob Faibussowitsch   PetscCall(DMGetMatType(dm, &mtype));
1679566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(*tdm, mtype));
168a0452a8eSToby Isaac   PetscFunctionReturn(0);
169a0452a8eSToby Isaac }
170a0452a8eSToby Isaac 
17101d9d024SToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm);
17201d9d024SToby Isaac 
173d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm)
174d71ae5a4SJacob Faibussowitsch {
175db4d5e8cSToby Isaac   DM_Forest  *forest = (DM_Forest *)dm->data;
176db4d5e8cSToby Isaac   const char *type;
177db4d5e8cSToby Isaac 
178db4d5e8cSToby Isaac   PetscFunctionBegin;
179db4d5e8cSToby Isaac   forest->refct++;
180db4d5e8cSToby Isaac   (*newdm)->data = forest;
1819566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetType((PetscObject)dm, &type));
1829566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, type));
1839566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Forest(*newdm));
184db4d5e8cSToby Isaac   PetscFunctionReturn(0);
185db4d5e8cSToby Isaac }
186db4d5e8cSToby Isaac 
187d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDestroy_Forest(DM dm)
188d71ae5a4SJacob Faibussowitsch {
189db4d5e8cSToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
190db4d5e8cSToby Isaac 
191db4d5e8cSToby Isaac   PetscFunctionBegin;
192db4d5e8cSToby Isaac   if (--forest->refct > 0) PetscFunctionReturn(0);
1939566063dSJacob Faibussowitsch   if (forest->destroy) PetscCall((*forest->destroy)(dm));
1949566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&forest->cellSF));
1959566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&forest->preCoarseToFine));
1969566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&forest->coarseToPreFine));
1979566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&forest->adaptLabel));
1989566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest->adaptStrategy));
1999566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&forest->base));
2009566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&forest->adapt));
2019566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest->topology));
2029566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest));
203db4d5e8cSToby Isaac   PetscFunctionReturn(0);
204db4d5e8cSToby Isaac }
205db4d5e8cSToby Isaac 
2069be51f97SToby Isaac /*@C
2079be51f97SToby Isaac   DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase.  The topology is a string (e.g.
20868d54884SBarry Smith   "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest during
2099be51f97SToby Isaac   DMSetUp().
2109be51f97SToby Isaac 
2119be51f97SToby Isaac   Logically collective on dm
2129be51f97SToby Isaac 
2139be51f97SToby Isaac   Input parameters:
2149be51f97SToby Isaac + dm - the forest
2159be51f97SToby Isaac - topology - the topology of the forest
2169be51f97SToby Isaac 
2179be51f97SToby Isaac   Level: intermediate
2189be51f97SToby Isaac 
219db781477SPatrick Sanan .seealso: `DMForestGetTopology()`, `DMForestSetBaseDM()`
2209be51f97SToby Isaac @*/
221d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology)
222d71ae5a4SJacob Faibussowitsch {
223db4d5e8cSToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
224db4d5e8cSToby Isaac 
225db4d5e8cSToby Isaac   PetscFunctionBegin;
226db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
22728b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the topology after setup");
2289566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest->topology));
2299566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy((const char *)topology, (char **)&forest->topology));
230db4d5e8cSToby Isaac   PetscFunctionReturn(0);
231db4d5e8cSToby Isaac }
232db4d5e8cSToby Isaac 
2339be51f97SToby Isaac /*@C
2349be51f97SToby Isaac   DMForestGetTopology - Get a string describing the topology of a DMForest.
2359be51f97SToby Isaac 
2369be51f97SToby Isaac   Not collective
2379be51f97SToby Isaac 
2389be51f97SToby Isaac   Input parameter:
2399be51f97SToby Isaac . dm - the forest
2409be51f97SToby Isaac 
2419be51f97SToby Isaac   Output parameter:
2429be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell')
2439be51f97SToby Isaac 
2449be51f97SToby Isaac   Level: intermediate
2459be51f97SToby Isaac 
246db781477SPatrick Sanan .seealso: `DMForestSetTopology()`
2479be51f97SToby Isaac @*/
248d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology)
249d71ae5a4SJacob Faibussowitsch {
250dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
251dd8e54a2SToby Isaac 
252dd8e54a2SToby Isaac   PetscFunctionBegin;
253dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
254dd8e54a2SToby Isaac   PetscValidPointer(topology, 2);
255dd8e54a2SToby Isaac   *topology = forest->topology;
256dd8e54a2SToby Isaac   PetscFunctionReturn(0);
257dd8e54a2SToby Isaac }
258dd8e54a2SToby Isaac 
2599be51f97SToby Isaac /*@
2609be51f97SToby Isaac   DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest.  The
2619be51f97SToby Isaac   forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its
262765b024eSBarry Smith   base.  In general, two forest must share a base to be comparable, to do things like construct interpolators.
2639be51f97SToby Isaac 
2649be51f97SToby Isaac   Logically collective on dm
2659be51f97SToby Isaac 
2669be51f97SToby Isaac   Input Parameters:
2679be51f97SToby Isaac + dm - the forest
2689be51f97SToby Isaac - base - the base DM of the forest
2699be51f97SToby Isaac 
27095452b02SPatrick Sanan   Notes:
27195452b02SPatrick Sanan     Currently the base DM must be a DMPLEX
272765b024eSBarry Smith 
2739be51f97SToby Isaac   Level: intermediate
2749be51f97SToby Isaac 
275db781477SPatrick Sanan .seealso: `DMForestGetBaseDM()`
2769be51f97SToby Isaac @*/
277d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetBaseDM(DM dm, DM base)
278d71ae5a4SJacob Faibussowitsch {
279dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
280dd8e54a2SToby Isaac   PetscInt   dim, dimEmbed;
281dd8e54a2SToby Isaac 
282dd8e54a2SToby Isaac   PetscFunctionBegin;
283dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
28428b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the base after setup");
2859566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)base));
2869566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&forest->base));
287dd8e54a2SToby Isaac   forest->base = base;
288a0452a8eSToby Isaac   if (base) {
2894fb89dddSMatthew G. Knepley     const PetscReal *maxCell, *Lstart, *L;
29028dfcf7cSStefano Zampini 
291a0452a8eSToby Isaac     PetscValidHeaderSpecific(base, DM_CLASSID, 2);
2929566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(base, &dim));
2939566063dSJacob Faibussowitsch     PetscCall(DMSetDimension(dm, dim));
2949566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDim(base, &dimEmbed));
2959566063dSJacob Faibussowitsch     PetscCall(DMSetCoordinateDim(dm, dimEmbed));
2964fb89dddSMatthew G. Knepley     PetscCall(DMGetPeriodicity(base, &maxCell, &Lstart, &L));
2974fb89dddSMatthew G. Knepley     PetscCall(DMSetPeriodicity(dm, maxCell, Lstart, L));
2984fb89dddSMatthew G. Knepley   } else PetscCall(DMSetPeriodicity(dm, NULL, NULL, NULL));
299dd8e54a2SToby Isaac   PetscFunctionReturn(0);
300dd8e54a2SToby Isaac }
301dd8e54a2SToby Isaac 
3029be51f97SToby Isaac /*@
3039be51f97SToby Isaac   DMForestGetBaseDM - Get the base DM of a DMForest forest.  The forest will be hierarchically refined from the base,
30468d54884SBarry Smith   and all refinements/coarsenings of the forest will share its base.  In general, two forest must share a base to be
3059be51f97SToby Isaac   comparable, to do things like construct interpolators.
3069be51f97SToby Isaac 
3079be51f97SToby Isaac   Not collective
3089be51f97SToby Isaac 
3099be51f97SToby Isaac   Input Parameter:
3109be51f97SToby Isaac . dm - the forest
3119be51f97SToby Isaac 
3129be51f97SToby Isaac   Output Parameter:
3139be51f97SToby Isaac . base - the base DM of the forest
3149be51f97SToby Isaac 
315367003a6SStefano Zampini   Notes:
316367003a6SStefano Zampini     After DMSetUp(), the base DM will be redundantly distributed across MPI processes
317367003a6SStefano Zampini 
3189be51f97SToby Isaac   Level: intermediate
3199be51f97SToby Isaac 
320db781477SPatrick Sanan .seealso: `DMForestSetBaseDM()`
3219be51f97SToby Isaac @*/
322d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetBaseDM(DM dm, DM *base)
323d71ae5a4SJacob Faibussowitsch {
324dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
325dd8e54a2SToby Isaac 
326dd8e54a2SToby Isaac   PetscFunctionBegin;
327dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
328dd8e54a2SToby Isaac   PetscValidPointer(base, 2);
329dd8e54a2SToby Isaac   *base = forest->base;
330dd8e54a2SToby Isaac   PetscFunctionReturn(0);
331dd8e54a2SToby Isaac }
332dd8e54a2SToby Isaac 
333d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *), void *ctx)
334d71ae5a4SJacob Faibussowitsch {
335cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
336cf38a08cSToby Isaac 
337cf38a08cSToby Isaac   PetscFunctionBegin;
338cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
339cf38a08cSToby Isaac   forest->mapcoordinates    = func;
340cf38a08cSToby Isaac   forest->mapcoordinatesctx = ctx;
341cf38a08cSToby Isaac   PetscFunctionReturn(0);
342cf38a08cSToby Isaac }
343cf38a08cSToby Isaac 
344d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *), void *ctx)
345d71ae5a4SJacob Faibussowitsch {
346cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
347cf38a08cSToby Isaac 
348cf38a08cSToby Isaac   PetscFunctionBegin;
349cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
350cf38a08cSToby Isaac   if (func) *func = forest->mapcoordinates;
351cf38a08cSToby Isaac   if (ctx) *((void **)ctx) = forest->mapcoordinatesctx;
352cf38a08cSToby Isaac   PetscFunctionReturn(0);
353cf38a08cSToby Isaac }
354cf38a08cSToby Isaac 
3559be51f97SToby Isaac /*@
3569be51f97SToby Isaac   DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be
3579be51f97SToby Isaac   adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp().  Usually not needed
3589be51f97SToby Isaac   by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this
3599be51f97SToby Isaac   routine.
3609be51f97SToby Isaac 
361dffe73a3SToby Isaac   Note that this can be called after setup with adapt = NULL, which will clear all internal data related to the
362dffe73a3SToby Isaac   adaptivity forest from dm.  This way, repeatedly adapting does not leave stale DM objects in memory.
363dffe73a3SToby Isaac 
3649be51f97SToby Isaac   Logically collective on dm
3659be51f97SToby Isaac 
366d8d19677SJose E. Roman   Input Parameters:
3679be51f97SToby Isaac + dm - the new forest, which will be constructed from adapt
3689be51f97SToby Isaac - adapt - the old forest
3699be51f97SToby Isaac 
3709be51f97SToby Isaac   Level: intermediate
3719be51f97SToby Isaac 
372db781477SPatrick Sanan .seealso: `DMForestGetAdaptivityForest()`, `DMForestSetAdaptivityPurpose()`
3739be51f97SToby Isaac @*/
374d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityForest(DM dm, DM adapt)
375d71ae5a4SJacob Faibussowitsch {
376dffe73a3SToby Isaac   DM_Forest *forest, *adaptForest, *oldAdaptForest;
377dffe73a3SToby Isaac   DM         oldAdapt;
378456cc5b7SMatthew G. Knepley   PetscBool  isForest;
379dd8e54a2SToby Isaac 
380dd8e54a2SToby Isaac   PetscFunctionBegin;
381dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3821fd50544SStefano Zampini   if (adapt) PetscValidHeaderSpecific(adapt, DM_CLASSID, 2);
3839566063dSJacob Faibussowitsch   PetscCall(DMIsForest(dm, &isForest));
384456cc5b7SMatthew G. Knepley   if (!isForest) PetscFunctionReturn(0);
3851dca8a05SBarry Smith   PetscCheck(adapt == NULL || !dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adaptation forest after setup");
386ba936b91SToby Isaac   forest = (DM_Forest *)dm->data;
3879566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivityForest(dm, &oldAdapt));
388193eb951SToby Isaac   adaptForest    = (DM_Forest *)(adapt ? adapt->data : NULL);
389193eb951SToby Isaac   oldAdaptForest = (DM_Forest *)(oldAdapt ? oldAdapt->data : NULL);
390dffe73a3SToby Isaac   if (adaptForest != oldAdaptForest) {
3919566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&forest->preCoarseToFine));
3929566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&forest->coarseToPreFine));
3939566063dSJacob Faibussowitsch     if (forest->clearadaptivityforest) PetscCall((*forest->clearadaptivityforest)(dm));
394dffe73a3SToby Isaac   }
39526d9498aSToby Isaac   switch (forest->adaptPurpose) {
396cd3c525cSToby Isaac   case DM_ADAPT_DETERMINE:
3979566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)adapt));
3989566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&(forest->adapt)));
399ba936b91SToby Isaac     forest->adapt = adapt;
40026d9498aSToby Isaac     break;
401d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_REFINE:
402d71ae5a4SJacob Faibussowitsch     PetscCall(DMSetCoarseDM(dm, adapt));
403d71ae5a4SJacob Faibussowitsch     break;
404a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
405d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_COARSEN_LAST:
406d71ae5a4SJacob Faibussowitsch     PetscCall(DMSetFineDM(dm, adapt));
407d71ae5a4SJacob Faibussowitsch     break;
408d71ae5a4SJacob Faibussowitsch   default:
409d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid adaptivity purpose");
41026d9498aSToby Isaac   }
411dd8e54a2SToby Isaac   PetscFunctionReturn(0);
412dd8e54a2SToby Isaac }
413dd8e54a2SToby Isaac 
4149be51f97SToby Isaac /*@
4159be51f97SToby Isaac   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.
4169be51f97SToby Isaac 
4179be51f97SToby Isaac   Not collective
4189be51f97SToby Isaac 
4199be51f97SToby Isaac   Input Parameter:
4209be51f97SToby Isaac . dm - the forest
4219be51f97SToby Isaac 
4229be51f97SToby Isaac   Output Parameter:
4239be51f97SToby Isaac . adapt - the forest from which dm is/was adapted
4249be51f97SToby Isaac 
4259be51f97SToby Isaac   Level: intermediate
4269be51f97SToby Isaac 
427db781477SPatrick Sanan .seealso: `DMForestSetAdaptivityForest()`, `DMForestSetAdaptivityPurpose()`
4289be51f97SToby Isaac @*/
429d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
430d71ae5a4SJacob Faibussowitsch {
431ba936b91SToby Isaac   DM_Forest *forest;
432dd8e54a2SToby Isaac 
433dd8e54a2SToby Isaac   PetscFunctionBegin;
434dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
435ba936b91SToby Isaac   forest = (DM_Forest *)dm->data;
43626d9498aSToby Isaac   switch (forest->adaptPurpose) {
437d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_DETERMINE:
438d71ae5a4SJacob Faibussowitsch     *adapt = forest->adapt;
439d71ae5a4SJacob Faibussowitsch     break;
440d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_REFINE:
441d71ae5a4SJacob Faibussowitsch     PetscCall(DMGetCoarseDM(dm, adapt));
442d71ae5a4SJacob Faibussowitsch     break;
443a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
444d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_COARSEN_LAST:
445d71ae5a4SJacob Faibussowitsch     PetscCall(DMGetFineDM(dm, adapt));
446d71ae5a4SJacob Faibussowitsch     break;
447d71ae5a4SJacob Faibussowitsch   default:
448d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid adaptivity purpose");
44926d9498aSToby Isaac   }
45026d9498aSToby Isaac   PetscFunctionReturn(0);
45126d9498aSToby Isaac }
45226d9498aSToby Isaac 
4539be51f97SToby Isaac /*@
4549be51f97SToby Isaac   DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its
455a1b0c543SToby Isaac   source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening
456cd3c525cSToby Isaac   (DM_ADAPT_COARSEN), or undefined (DM_ADAPT_DETERMINE).  This only matters for the purposes of reference counting:
4579be51f97SToby Isaac   during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse
4589be51f97SToby Isaac   relationship (see DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does
4599be51f97SToby Isaac   not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can
4609be51f97SToby Isaac   cause a memory leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
4619be51f97SToby Isaac 
4629be51f97SToby Isaac   Logically collective on dm
4639be51f97SToby Isaac 
4649be51f97SToby Isaac   Input Parameters:
4659be51f97SToby Isaac + dm - the forest
466bf2d5fbbSStefano Zampini - purpose - the adaptivity purpose
4679be51f97SToby Isaac 
4689be51f97SToby Isaac   Level: advanced
4699be51f97SToby Isaac 
470db781477SPatrick Sanan .seealso: `DMForestTemplate()`, `DMForestSetAdaptivityForest()`, `DMForestGetAdaptivityForest()`, `DMAdaptFlag`
4719be51f97SToby Isaac @*/
472d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose)
473d71ae5a4SJacob Faibussowitsch {
47426d9498aSToby Isaac   DM_Forest *forest;
47526d9498aSToby Isaac 
47626d9498aSToby Isaac   PetscFunctionBegin;
47726d9498aSToby Isaac   forest = (DM_Forest *)dm->data;
47828b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adaptation forest after setup");
47926d9498aSToby Isaac   if (purpose != forest->adaptPurpose) {
48026d9498aSToby Isaac     DM adapt;
48126d9498aSToby Isaac 
4829566063dSJacob Faibussowitsch     PetscCall(DMForestGetAdaptivityForest(dm, &adapt));
4839566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)adapt));
4849566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, NULL));
485f885a11aSToby Isaac 
48626d9498aSToby Isaac     forest->adaptPurpose = purpose;
487f885a11aSToby Isaac 
4889566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, adapt));
4899566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&adapt));
49026d9498aSToby Isaac   }
491dd8e54a2SToby Isaac   PetscFunctionReturn(0);
492dd8e54a2SToby Isaac }
493dd8e54a2SToby Isaac 
49456c0450aSToby Isaac /*@
49556c0450aSToby Isaac   DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with
496bf2d5fbbSStefano Zampini   DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening (DM_ADAPT_COARSEN),
497bf2d5fbbSStefano Zampini   coarsening only the last level (DM_ADAPT_COARSEN_LAST) or undefined (DM_ADAPT_DETERMINE).
498bf2d5fbbSStefano Zampini   This only matters for the purposes of reference counting: during DMDestroy(), cyclic
49956c0450aSToby Isaac   references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see
50056c0450aSToby Isaac   DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does not maintain a
50156c0450aSToby Isaac   reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory
50256c0450aSToby Isaac   leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
50356c0450aSToby Isaac 
50456c0450aSToby Isaac   Not collective
50556c0450aSToby Isaac 
50656c0450aSToby Isaac   Input Parameter:
50756c0450aSToby Isaac . dm - the forest
50856c0450aSToby Isaac 
50956c0450aSToby Isaac   Output Parameter:
510bf2d5fbbSStefano Zampini . purpose - the adaptivity purpose
51156c0450aSToby Isaac 
51256c0450aSToby Isaac   Level: advanced
51356c0450aSToby Isaac 
514db781477SPatrick Sanan .seealso: `DMForestTemplate()`, `DMForestSetAdaptivityForest()`, `DMForestGetAdaptivityForest()`, `DMAdaptFlag`
51556c0450aSToby Isaac @*/
516d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose)
517d71ae5a4SJacob Faibussowitsch {
51856c0450aSToby Isaac   DM_Forest *forest;
51956c0450aSToby Isaac 
52056c0450aSToby Isaac   PetscFunctionBegin;
52156c0450aSToby Isaac   forest   = (DM_Forest *)dm->data;
52256c0450aSToby Isaac   *purpose = forest->adaptPurpose;
52356c0450aSToby Isaac   PetscFunctionReturn(0);
52456c0450aSToby Isaac }
52556c0450aSToby Isaac 
5269be51f97SToby Isaac /*@
5279be51f97SToby Isaac   DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
5289be51f97SToby Isaac   cell adjacency (for the purposes of partitioning and overlap).
5299be51f97SToby Isaac 
5309be51f97SToby Isaac   Logically collective on dm
5319be51f97SToby Isaac 
5329be51f97SToby Isaac   Input Parameters:
5339be51f97SToby Isaac + dm - the forest
5349be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency)
5359be51f97SToby Isaac 
5369be51f97SToby Isaac   Level: intermediate
5379be51f97SToby Isaac 
538db781477SPatrick Sanan .seealso: `DMForestGetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`, `DMForestSetPartitionOverlap()`
5399be51f97SToby Isaac @*/
540d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
541d71ae5a4SJacob Faibussowitsch {
542dd8e54a2SToby Isaac   PetscInt   dim;
543dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
544dd8e54a2SToby Isaac 
545dd8e54a2SToby Isaac   PetscFunctionBegin;
546dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
54728b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adjacency dimension after setup");
54863a3b9bcSJacob Faibussowitsch   PetscCheck(adjDim >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "adjacency dim cannot be < 0: %" PetscInt_FMT, adjDim);
5499566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
55063a3b9bcSJacob Faibussowitsch   PetscCheck(adjDim <= dim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "adjacency dim cannot be > %" PetscInt_FMT ": %" PetscInt_FMT, dim, adjDim);
551dd8e54a2SToby Isaac   forest->adjDim = adjDim;
552dd8e54a2SToby Isaac   PetscFunctionReturn(0);
553dd8e54a2SToby Isaac }
554dd8e54a2SToby Isaac 
5559be51f97SToby Isaac /*@
5569be51f97SToby Isaac   DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that,
5579be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
5589be51f97SToby Isaac 
5599be51f97SToby Isaac   Logically collective on dm
5609be51f97SToby Isaac 
5619be51f97SToby Isaac   Input Parameters:
5629be51f97SToby Isaac + dm - the forest
5639be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
5649be51f97SToby Isaac 
5659be51f97SToby Isaac   Level: intermediate
5669be51f97SToby Isaac 
567db781477SPatrick Sanan .seealso: `DMForestGetAdjacencyCodimension()`, `DMForestSetAdjacencyDimension()`
5689be51f97SToby Isaac @*/
569d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
570d71ae5a4SJacob Faibussowitsch {
571dd8e54a2SToby Isaac   PetscInt dim;
572dd8e54a2SToby Isaac 
573dd8e54a2SToby Isaac   PetscFunctionBegin;
574dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5759566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
5769566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdjacencyDimension(dm, dim - adjCodim));
577dd8e54a2SToby Isaac   PetscFunctionReturn(0);
578dd8e54a2SToby Isaac }
579dd8e54a2SToby Isaac 
5809be51f97SToby Isaac /*@
5819be51f97SToby Isaac   DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
5829be51f97SToby Isaac   purposes of partitioning and overlap).
5839be51f97SToby Isaac 
5849be51f97SToby Isaac   Not collective
5859be51f97SToby Isaac 
5869be51f97SToby Isaac   Input Parameter:
5879be51f97SToby Isaac . dm - the forest
5889be51f97SToby Isaac 
5899be51f97SToby Isaac   Output Parameter:
5909be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency)
5919be51f97SToby Isaac 
5929be51f97SToby Isaac   Level: intermediate
5939be51f97SToby Isaac 
594db781477SPatrick Sanan .seealso: `DMForestSetAdjacencyDimension()`, `DMForestGetAdjacencyCodimension()`, `DMForestSetPartitionOverlap()`
5959be51f97SToby Isaac @*/
596d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
597d71ae5a4SJacob Faibussowitsch {
598dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
599dd8e54a2SToby Isaac 
600dd8e54a2SToby Isaac   PetscFunctionBegin;
601dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
602dd8e54a2SToby Isaac   PetscValidIntPointer(adjDim, 2);
603dd8e54a2SToby Isaac   *adjDim = forest->adjDim;
604dd8e54a2SToby Isaac   PetscFunctionReturn(0);
605dd8e54a2SToby Isaac }
606dd8e54a2SToby Isaac 
6079be51f97SToby Isaac /*@
6089be51f97SToby Isaac   DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that,
6099be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
6109be51f97SToby Isaac 
6119be51f97SToby Isaac   Not collective
6129be51f97SToby Isaac 
6139be51f97SToby Isaac   Input Parameter:
6149be51f97SToby Isaac . dm - the forest
6159be51f97SToby Isaac 
6169be51f97SToby Isaac   Output Parameter:
6179be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
6189be51f97SToby Isaac 
6199be51f97SToby Isaac   Level: intermediate
6209be51f97SToby Isaac 
621db781477SPatrick Sanan .seealso: `DMForestSetAdjacencyCodimension()`, `DMForestGetAdjacencyDimension()`
6229be51f97SToby Isaac @*/
623d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
624d71ae5a4SJacob Faibussowitsch {
625dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
626dd8e54a2SToby Isaac   PetscInt   dim;
627dd8e54a2SToby Isaac 
628dd8e54a2SToby Isaac   PetscFunctionBegin;
629dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
630dd8e54a2SToby Isaac   PetscValidIntPointer(adjCodim, 2);
6319566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
632dd8e54a2SToby Isaac   *adjCodim = dim - forest->adjDim;
633dd8e54a2SToby Isaac   PetscFunctionReturn(0);
634dd8e54a2SToby Isaac }
635dd8e54a2SToby Isaac 
6369be51f97SToby Isaac /*@
6379be51f97SToby Isaac   DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel
6389be51f97SToby Isaac   partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding
6399be51f97SToby Isaac   adjacent cells
6409be51f97SToby Isaac 
6419be51f97SToby Isaac   Logically collective on dm
6429be51f97SToby Isaac 
6439be51f97SToby Isaac   Input Parameters:
6449be51f97SToby Isaac + dm - the forest
6459be51f97SToby Isaac - overlap - default 0
6469be51f97SToby Isaac 
6479be51f97SToby Isaac   Level: intermediate
6489be51f97SToby Isaac 
649db781477SPatrick Sanan .seealso: `DMForestGetPartitionOverlap()`, `DMForestSetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`
6509be51f97SToby Isaac @*/
651d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap)
652d71ae5a4SJacob Faibussowitsch {
653dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
654dd8e54a2SToby Isaac 
655dd8e54a2SToby Isaac   PetscFunctionBegin;
656dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
65728b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the overlap after setup");
65863a3b9bcSJacob Faibussowitsch   PetscCheck(overlap >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "overlap cannot be < 0: %" PetscInt_FMT, overlap);
659dd8e54a2SToby Isaac   forest->overlap = overlap;
660dd8e54a2SToby Isaac   PetscFunctionReturn(0);
661dd8e54a2SToby Isaac }
662dd8e54a2SToby Isaac 
6639be51f97SToby Isaac /*@
6649be51f97SToby Isaac   DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values
6659be51f97SToby Isaac   > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells
6669be51f97SToby Isaac 
6679be51f97SToby Isaac   Not collective
6689be51f97SToby Isaac 
6699be51f97SToby Isaac   Input Parameter:
6709be51f97SToby Isaac . dm - the forest
6719be51f97SToby Isaac 
6729be51f97SToby Isaac   Output Parameter:
6739be51f97SToby Isaac . overlap - default 0
6749be51f97SToby Isaac 
6759be51f97SToby Isaac   Level: intermediate
6769be51f97SToby Isaac 
677db781477SPatrick Sanan .seealso: `DMForestGetPartitionOverlap()`, `DMForestSetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`
6789be51f97SToby Isaac @*/
679d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap)
680d71ae5a4SJacob Faibussowitsch {
681dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
682dd8e54a2SToby Isaac 
683dd8e54a2SToby Isaac   PetscFunctionBegin;
684dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
685dd8e54a2SToby Isaac   PetscValidIntPointer(overlap, 2);
686dd8e54a2SToby Isaac   *overlap = forest->overlap;
687dd8e54a2SToby Isaac   PetscFunctionReturn(0);
688dd8e54a2SToby Isaac }
689dd8e54a2SToby Isaac 
6909be51f97SToby Isaac /*@
6919be51f97SToby Isaac   DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base
6929be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest
6939be51f97SToby Isaac   (see DMForestGetAdaptivityForest()) this limits the amount of coarsening.
6949be51f97SToby Isaac 
6959be51f97SToby Isaac   Logically collective on dm
6969be51f97SToby Isaac 
6979be51f97SToby Isaac   Input Parameters:
6989be51f97SToby Isaac + dm - the forest
6999be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7009be51f97SToby Isaac 
7019be51f97SToby Isaac   Level: intermediate
7029be51f97SToby Isaac 
703db781477SPatrick Sanan .seealso: `DMForestGetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestSetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()`
7049be51f97SToby Isaac @*/
705d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement)
706d71ae5a4SJacob Faibussowitsch {
707dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
708dd8e54a2SToby Isaac 
709dd8e54a2SToby Isaac   PetscFunctionBegin;
710dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
71128b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the minimum refinement after setup");
712dd8e54a2SToby Isaac   forest->minRefinement = minRefinement;
713dd8e54a2SToby Isaac   PetscFunctionReturn(0);
714dd8e54a2SToby Isaac }
715dd8e54a2SToby Isaac 
7169be51f97SToby Isaac /*@
7179be51f97SToby Isaac   DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see
7189be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest (see
7199be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of coarsening.
7209be51f97SToby Isaac 
7219be51f97SToby Isaac   Not collective
7229be51f97SToby Isaac 
7239be51f97SToby Isaac   Input Parameter:
7249be51f97SToby Isaac . dm - the forest
7259be51f97SToby Isaac 
7269be51f97SToby Isaac   Output Parameter:
7279be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7289be51f97SToby Isaac 
7299be51f97SToby Isaac   Level: intermediate
7309be51f97SToby Isaac 
731db781477SPatrick Sanan .seealso: `DMForestSetMinimumRefinement()`, `DMForestGetMaximumRefinement()`, `DMForestGetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()`
7329be51f97SToby Isaac @*/
733d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement)
734d71ae5a4SJacob Faibussowitsch {
735dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
736dd8e54a2SToby Isaac 
737dd8e54a2SToby Isaac   PetscFunctionBegin;
738dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
739dd8e54a2SToby Isaac   PetscValidIntPointer(minRefinement, 2);
740dd8e54a2SToby Isaac   *minRefinement = forest->minRefinement;
741dd8e54a2SToby Isaac   PetscFunctionReturn(0);
742dd8e54a2SToby Isaac }
743dd8e54a2SToby Isaac 
7449be51f97SToby Isaac /*@
7459be51f97SToby Isaac   DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base
7469be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.
7479be51f97SToby Isaac 
7489be51f97SToby Isaac   Logically collective on dm
7499be51f97SToby Isaac 
7509be51f97SToby Isaac   Input Parameters:
7519be51f97SToby Isaac + dm - the forest
7529be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7539be51f97SToby Isaac 
7549be51f97SToby Isaac   Level: intermediate
7559be51f97SToby Isaac 
756db781477SPatrick Sanan .seealso: `DMForestSetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestGetBaseDM()`
7579be51f97SToby Isaac @*/
758d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement)
759d71ae5a4SJacob Faibussowitsch {
76056ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
76156ba9f64SToby Isaac 
76256ba9f64SToby Isaac   PetscFunctionBegin;
76356ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
76428b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the initial refinement after setup");
76556ba9f64SToby Isaac   forest->initRefinement = initRefinement;
76656ba9f64SToby Isaac   PetscFunctionReturn(0);
76756ba9f64SToby Isaac }
76856ba9f64SToby Isaac 
7699be51f97SToby Isaac /*@
7709be51f97SToby Isaac   DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see
7719be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.
7729be51f97SToby Isaac 
7739be51f97SToby Isaac   Not collective
7749be51f97SToby Isaac 
7759be51f97SToby Isaac   Input Parameter:
7769be51f97SToby Isaac . dm - the forest
7779be51f97SToby Isaac 
77801d2d390SJose E. Roman   Output Parameter:
779bf2d5fbbSStefano Zampini . initRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7809be51f97SToby Isaac 
7819be51f97SToby Isaac   Level: intermediate
7829be51f97SToby Isaac 
783db781477SPatrick Sanan .seealso: `DMForestSetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestGetBaseDM()`
7849be51f97SToby Isaac @*/
785d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement)
786d71ae5a4SJacob Faibussowitsch {
78756ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
78856ba9f64SToby Isaac 
78956ba9f64SToby Isaac   PetscFunctionBegin;
79056ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
79156ba9f64SToby Isaac   PetscValidIntPointer(initRefinement, 2);
79256ba9f64SToby Isaac   *initRefinement = forest->initRefinement;
79356ba9f64SToby Isaac   PetscFunctionReturn(0);
79456ba9f64SToby Isaac }
79556ba9f64SToby Isaac 
7969be51f97SToby Isaac /*@
7979be51f97SToby Isaac   DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base
7989be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest
7999be51f97SToby Isaac   (see DMForestGetAdaptivityForest()), this limits the amount of refinement.
8009be51f97SToby Isaac 
8019be51f97SToby Isaac   Logically collective on dm
8029be51f97SToby Isaac 
8039be51f97SToby Isaac   Input Parameters:
8049be51f97SToby Isaac + dm - the forest
8059be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8069be51f97SToby Isaac 
8079be51f97SToby Isaac   Level: intermediate
8089be51f97SToby Isaac 
809db781477SPatrick Sanan .seealso: `DMForestGetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestSetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityDM()`
8109be51f97SToby Isaac @*/
811d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement)
812d71ae5a4SJacob Faibussowitsch {
813dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
814dd8e54a2SToby Isaac 
815dd8e54a2SToby Isaac   PetscFunctionBegin;
816dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
81728b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the maximum refinement after setup");
818c7eeac06SToby Isaac   forest->maxRefinement = maxRefinement;
819dd8e54a2SToby Isaac   PetscFunctionReturn(0);
820dd8e54a2SToby Isaac }
821dd8e54a2SToby Isaac 
8229be51f97SToby Isaac /*@
8239be51f97SToby Isaac   DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see
8249be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest (see
8259be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of refinement.
8269be51f97SToby Isaac 
8279be51f97SToby Isaac   Not collective
8289be51f97SToby Isaac 
8299be51f97SToby Isaac   Input Parameter:
8309be51f97SToby Isaac . dm - the forest
8319be51f97SToby Isaac 
8329be51f97SToby Isaac   Output Parameter:
8339be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8349be51f97SToby Isaac 
8359be51f97SToby Isaac   Level: intermediate
8369be51f97SToby Isaac 
837db781477SPatrick Sanan .seealso: `DMForestSetMaximumRefinement()`, `DMForestGetMinimumRefinement()`, `DMForestGetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()`
8389be51f97SToby Isaac @*/
839d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement)
840d71ae5a4SJacob Faibussowitsch {
841dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
842dd8e54a2SToby Isaac 
843dd8e54a2SToby Isaac   PetscFunctionBegin;
844dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
845c7eeac06SToby Isaac   PetscValidIntPointer(maxRefinement, 2);
846c7eeac06SToby Isaac   *maxRefinement = forest->maxRefinement;
847dd8e54a2SToby Isaac   PetscFunctionReturn(0);
848dd8e54a2SToby Isaac }
849c7eeac06SToby Isaac 
8509be51f97SToby Isaac /*@C
8519be51f97SToby Isaac   DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes.
8529be51f97SToby Isaac   Subtypes of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree
8539be51f97SToby Isaac   for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to
8549be51f97SToby Isaac   specify refinement/coarsening.
8559be51f97SToby Isaac 
8569be51f97SToby Isaac   Logically collective on dm
8579be51f97SToby Isaac 
8589be51f97SToby Isaac   Input Parameters:
8599be51f97SToby Isaac + dm - the forest
8609be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL
8619be51f97SToby Isaac 
8629be51f97SToby Isaac   Level: advanced
8639be51f97SToby Isaac 
864db781477SPatrick Sanan .seealso: `DMForestGetAdaptivityStrategy()`
8659be51f97SToby Isaac @*/
866d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy)
867d71ae5a4SJacob Faibussowitsch {
868c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
869c7eeac06SToby Isaac 
870c7eeac06SToby Isaac   PetscFunctionBegin;
871c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8729566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest->adaptStrategy));
8739566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy((const char *)adaptStrategy, (char **)&forest->adaptStrategy));
874c7eeac06SToby Isaac   PetscFunctionReturn(0);
875c7eeac06SToby Isaac }
876c7eeac06SToby Isaac 
8779be51f97SToby Isaac /*@C
8789be51f97SToby Isaac   DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes.  Subtypes
8799be51f97SToby Isaac   of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all
8809be51f97SToby Isaac   processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only
8819be51f97SToby Isaac   one process needs to specify refinement/coarsening.
8829be51f97SToby Isaac 
8839be51f97SToby Isaac   Not collective
8849be51f97SToby Isaac 
8859be51f97SToby Isaac   Input Parameter:
8869be51f97SToby Isaac . dm - the forest
8879be51f97SToby Isaac 
8889be51f97SToby Isaac   Output Parameter:
8899be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL)
8909be51f97SToby Isaac 
8919be51f97SToby Isaac   Level: advanced
8929be51f97SToby Isaac 
893db781477SPatrick Sanan .seealso: `DMForestSetAdaptivityStrategy()`
8949be51f97SToby Isaac @*/
895d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
896d71ae5a4SJacob Faibussowitsch {
897c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
898c7eeac06SToby Isaac 
899c7eeac06SToby Isaac   PetscFunctionBegin;
900c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
901c7eeac06SToby Isaac   PetscValidPointer(adaptStrategy, 2);
902c7eeac06SToby Isaac   *adaptStrategy = forest->adaptStrategy;
903c7eeac06SToby Isaac   PetscFunctionReturn(0);
904c7eeac06SToby Isaac }
905c7eeac06SToby Isaac 
9062a133e43SToby Isaac /*@
9072a133e43SToby Isaac   DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning,
9082a133e43SToby Isaac   etc.) was successful.  PETSC_FALSE indicates that the post-adaptation forest is the same as the pre-adpatation
9092a133e43SToby Isaac   forest.  A requested adaptation may have been unsuccessful if, for example, the requested refinement would have
9102a133e43SToby Isaac   exceeded the maximum refinement level.
9112a133e43SToby Isaac 
9122a133e43SToby Isaac   Collective on dm
9132a133e43SToby Isaac 
9142a133e43SToby Isaac   Input Parameter:
9152a133e43SToby Isaac 
9162a133e43SToby Isaac . dm - the post-adaptation forest
9172a133e43SToby Isaac 
9182a133e43SToby Isaac   Output Parameter:
9192a133e43SToby Isaac 
9202a133e43SToby Isaac . success - PETSC_TRUE if the post-adaptation forest is different from the pre-adaptation forest.
9212a133e43SToby Isaac 
9222a133e43SToby Isaac   Level: intermediate
9232a133e43SToby Isaac 
9242a133e43SToby Isaac .see
9252a133e43SToby Isaac @*/
926d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success)
927d71ae5a4SJacob Faibussowitsch {
9282a133e43SToby Isaac   DM_Forest *forest;
9292a133e43SToby Isaac 
9302a133e43SToby Isaac   PetscFunctionBegin;
9312a133e43SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
93228b400f6SJacob Faibussowitsch   PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() has not been called yet.");
9332a133e43SToby Isaac   forest = (DM_Forest *)dm->data;
9349566063dSJacob Faibussowitsch   PetscCall((forest->getadaptivitysuccess)(dm, success));
9352a133e43SToby Isaac   PetscFunctionReturn(0);
9362a133e43SToby Isaac }
9372a133e43SToby Isaac 
938bf9b5d84SToby Isaac /*@
939bf9b5d84SToby Isaac   DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed
940bf9b5d84SToby Isaac   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().
941bf9b5d84SToby Isaac 
942bf9b5d84SToby Isaac   Logically collective on dm
943bf9b5d84SToby Isaac 
944bf9b5d84SToby Isaac   Input Parameters:
945bf9b5d84SToby Isaac + dm - the post-adaptation forest
946bf9b5d84SToby Isaac - computeSF - default PETSC_TRUE
947bf9b5d84SToby Isaac 
948bf9b5d84SToby Isaac   Level: advanced
949bf9b5d84SToby Isaac 
950db781477SPatrick Sanan .seealso: `DMForestGetComputeAdaptivitySF()`, `DMForestGetAdaptivitySF()`
951bf9b5d84SToby Isaac @*/
952d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
953d71ae5a4SJacob Faibussowitsch {
954bf9b5d84SToby Isaac   DM_Forest *forest;
955bf9b5d84SToby Isaac 
956bf9b5d84SToby Isaac   PetscFunctionBegin;
957bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
95828b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot compute adaptivity PetscSFs after setup is called");
959bf9b5d84SToby Isaac   forest                 = (DM_Forest *)dm->data;
960bf9b5d84SToby Isaac   forest->computeAdaptSF = computeSF;
961bf9b5d84SToby Isaac   PetscFunctionReturn(0);
962bf9b5d84SToby Isaac }
963bf9b5d84SToby Isaac 
964d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time)
965d71ae5a4SJacob Faibussowitsch {
96680b27e07SToby Isaac   DM_Forest *forest;
96780b27e07SToby Isaac 
96880b27e07SToby Isaac   PetscFunctionBegin;
96980b27e07SToby Isaac   PetscValidHeaderSpecific(dmIn, DM_CLASSID, 1);
97080b27e07SToby Isaac   PetscValidHeaderSpecific(vecIn, VEC_CLASSID, 2);
97180b27e07SToby Isaac   PetscValidHeaderSpecific(dmOut, DM_CLASSID, 3);
97280b27e07SToby Isaac   PetscValidHeaderSpecific(vecOut, VEC_CLASSID, 4);
97380b27e07SToby Isaac   forest = (DM_Forest *)dmIn->data;
97428b400f6SJacob Faibussowitsch   PetscCheck(forest->transfervec, PetscObjectComm((PetscObject)dmIn), PETSC_ERR_SUP, "DMForestTransferVec() not implemented");
9759566063dSJacob Faibussowitsch   PetscCall((forest->transfervec)(dmIn, vecIn, dmOut, vecOut, useBCs, time));
97680b27e07SToby Isaac   PetscFunctionReturn(0);
97780b27e07SToby Isaac }
97880b27e07SToby Isaac 
979d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestTransferVecFromBase(DM dm, Vec vecIn, Vec vecOut)
980d71ae5a4SJacob Faibussowitsch {
981ac34a06fSStefano Zampini   DM_Forest *forest;
982ac34a06fSStefano Zampini 
983ac34a06fSStefano Zampini   PetscFunctionBegin;
984ac34a06fSStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
985ac34a06fSStefano Zampini   PetscValidHeaderSpecific(vecIn, VEC_CLASSID, 2);
986ac34a06fSStefano Zampini   PetscValidHeaderSpecific(vecOut, VEC_CLASSID, 3);
987ac34a06fSStefano Zampini   forest = (DM_Forest *)dm->data;
98828b400f6SJacob Faibussowitsch   PetscCheck(forest->transfervecfrombase, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMForestTransferVecFromBase() not implemented");
9899566063dSJacob Faibussowitsch   PetscCall((forest->transfervecfrombase)(dm, vecIn, vecOut));
990ac34a06fSStefano Zampini   PetscFunctionReturn(0);
991ac34a06fSStefano Zampini }
992ac34a06fSStefano Zampini 
993bf9b5d84SToby Isaac /*@
994bf9b5d84SToby Isaac   DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the
995bf9b5d84SToby Isaac   pre-adaptation forest to the post-adaptiation forest.  After DMSetUp() is called, these transfer PetscSFs can be
996bf9b5d84SToby Isaac   accessed with DMForestGetAdaptivitySF().
997bf9b5d84SToby Isaac 
998bf9b5d84SToby Isaac   Not collective
999bf9b5d84SToby Isaac 
1000bf9b5d84SToby Isaac   Input Parameter:
1001bf9b5d84SToby Isaac . dm - the post-adaptation forest
1002bf9b5d84SToby Isaac 
1003bf9b5d84SToby Isaac   Output Parameter:
1004bf9b5d84SToby Isaac . computeSF - default PETSC_TRUE
1005bf9b5d84SToby Isaac 
1006bf9b5d84SToby Isaac   Level: advanced
1007bf9b5d84SToby Isaac 
1008db781477SPatrick Sanan .seealso: `DMForestSetComputeAdaptivitySF()`, `DMForestGetAdaptivitySF()`
1009bf9b5d84SToby Isaac @*/
1010d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1011d71ae5a4SJacob Faibussowitsch {
1012bf9b5d84SToby Isaac   DM_Forest *forest;
1013bf9b5d84SToby Isaac 
1014bf9b5d84SToby Isaac   PetscFunctionBegin;
1015bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1016bf9b5d84SToby Isaac   forest     = (DM_Forest *)dm->data;
1017bf9b5d84SToby Isaac   *computeSF = forest->computeAdaptSF;
1018bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1019bf9b5d84SToby Isaac }
1020bf9b5d84SToby Isaac 
1021bf9b5d84SToby Isaac /*@
1022bf9b5d84SToby Isaac   DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest.
1023bf9b5d84SToby Isaac   Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be
1024bf9b5d84SToby Isaac   some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa.  Therefore there are two
1025bf9b5d84SToby Isaac   PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates
1026bf9b5d84SToby Isaac   pre-adaptation fine cells to post-adaptation coarse cells.
1027bf9b5d84SToby Isaac 
1028bf9b5d84SToby Isaac   Not collective
1029bf9b5d84SToby Isaac 
1030bf9b5d84SToby Isaac   Input Parameter:
1031bf9b5d84SToby Isaac   dm - the post-adaptation forest
1032bf9b5d84SToby Isaac 
1033bf9b5d84SToby Isaac   Output Parameter:
10340f17b9e3SToby Isaac   preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
10350f17b9e3SToby Isaac   coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-
1036bf9b5d84SToby Isaac 
1037bf9b5d84SToby Isaac   Level: advanced
1038bf9b5d84SToby Isaac 
1039db781477SPatrick Sanan .seealso: `DMForestGetComputeAdaptivitySF()`, `DMForestSetComputeAdaptivitySF()`
1040bf9b5d84SToby Isaac @*/
1041d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1042d71ae5a4SJacob Faibussowitsch {
1043bf9b5d84SToby Isaac   DM_Forest *forest;
1044bf9b5d84SToby Isaac 
1045bf9b5d84SToby Isaac   PetscFunctionBegin;
1046bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10479566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
1048bf9b5d84SToby Isaac   forest = (DM_Forest *)dm->data;
1049f885a11aSToby Isaac   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1050f885a11aSToby Isaac   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
1051bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1052bf9b5d84SToby Isaac }
1053bf9b5d84SToby Isaac 
10549be51f97SToby Isaac /*@
10559be51f97SToby Isaac   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
10569be51f97SToby Isaac   indicate that the diameter of neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may
10579be51f97SToby Isaac   only support one particular choice of grading factor.
10589be51f97SToby Isaac 
10599be51f97SToby Isaac   Logically collective on dm
10609be51f97SToby Isaac 
10619be51f97SToby Isaac   Input Parameters:
10629be51f97SToby Isaac + dm - the forest
10639be51f97SToby Isaac - grade - the grading factor
10649be51f97SToby Isaac 
10659be51f97SToby Isaac   Level: advanced
10669be51f97SToby Isaac 
1067db781477SPatrick Sanan .seealso: `DMForestGetGradeFactor()`
10689be51f97SToby Isaac @*/
1069d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1070d71ae5a4SJacob Faibussowitsch {
1071c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1072c7eeac06SToby Isaac 
1073c7eeac06SToby Isaac   PetscFunctionBegin;
1074c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
107528b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the grade factor after setup");
1076c7eeac06SToby Isaac   forest->gradeFactor = grade;
1077c7eeac06SToby Isaac   PetscFunctionReturn(0);
1078c7eeac06SToby Isaac }
1079c7eeac06SToby Isaac 
10809be51f97SToby Isaac /*@
10819be51f97SToby Isaac   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
10829be51f97SToby Isaac   neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may only support one particular
10839be51f97SToby Isaac   choice of grading factor.
10849be51f97SToby Isaac 
10859be51f97SToby Isaac   Not collective
10869be51f97SToby Isaac 
10879be51f97SToby Isaac   Input Parameter:
10889be51f97SToby Isaac . dm - the forest
10899be51f97SToby Isaac 
10909be51f97SToby Isaac   Output Parameter:
10919be51f97SToby Isaac . grade - the grading factor
10929be51f97SToby Isaac 
10939be51f97SToby Isaac   Level: advanced
10949be51f97SToby Isaac 
1095db781477SPatrick Sanan .seealso: `DMForestSetGradeFactor()`
10969be51f97SToby Isaac @*/
1097d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1098d71ae5a4SJacob Faibussowitsch {
1099c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1100c7eeac06SToby Isaac 
1101c7eeac06SToby Isaac   PetscFunctionBegin;
1102c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1103c7eeac06SToby Isaac   PetscValidIntPointer(grade, 2);
1104c7eeac06SToby Isaac   *grade = forest->gradeFactor;
1105c7eeac06SToby Isaac   PetscFunctionReturn(0);
1106c7eeac06SToby Isaac }
1107c7eeac06SToby Isaac 
11089be51f97SToby Isaac /*@
11099be51f97SToby Isaac   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
11109be51f97SToby Isaac   the cell weight (see DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be
11119be51f97SToby Isaac   (cellWeight) * (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on
11129be51f97SToby Isaac   its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
11139be51f97SToby Isaac   computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11149be51f97SToby Isaac 
11159be51f97SToby Isaac   Logically collective on dm
11169be51f97SToby Isaac 
11179be51f97SToby Isaac   Input Parameters:
11189be51f97SToby Isaac + dm - the forest
11199be51f97SToby Isaac - weightsFactors - default 1.
11209be51f97SToby Isaac 
11219be51f97SToby Isaac   Level: advanced
11229be51f97SToby Isaac 
1123db781477SPatrick Sanan .seealso: `DMForestGetCellWeightFactor()`, `DMForestSetCellWeights()`
11249be51f97SToby Isaac @*/
1125d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
1126d71ae5a4SJacob Faibussowitsch {
1127c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1128c7eeac06SToby Isaac 
1129c7eeac06SToby Isaac   PetscFunctionBegin;
1130c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
113128b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the weights factor after setup");
1132c7eeac06SToby Isaac   forest->weightsFactor = weightsFactor;
1133c7eeac06SToby Isaac   PetscFunctionReturn(0);
1134c7eeac06SToby Isaac }
1135c7eeac06SToby Isaac 
11369be51f97SToby Isaac /*@
11379be51f97SToby Isaac   DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
11389be51f97SToby Isaac   DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be (cellWeight) *
11399be51f97SToby Isaac   (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on its level; a
11409be51f97SToby Isaac   factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation
11419be51f97SToby Isaac   associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11429be51f97SToby Isaac 
11439be51f97SToby Isaac   Not collective
11449be51f97SToby Isaac 
11459be51f97SToby Isaac   Input Parameter:
11469be51f97SToby Isaac . dm - the forest
11479be51f97SToby Isaac 
11489be51f97SToby Isaac   Output Parameter:
11499be51f97SToby Isaac . weightsFactors - default 1.
11509be51f97SToby Isaac 
11519be51f97SToby Isaac   Level: advanced
11529be51f97SToby Isaac 
1153db781477SPatrick Sanan .seealso: `DMForestSetCellWeightFactor()`, `DMForestSetCellWeights()`
11549be51f97SToby Isaac @*/
1155d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1156d71ae5a4SJacob Faibussowitsch {
1157c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1158c7eeac06SToby Isaac 
1159c7eeac06SToby Isaac   PetscFunctionBegin;
1160c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1161c7eeac06SToby Isaac   PetscValidRealPointer(weightsFactor, 2);
1162c7eeac06SToby Isaac   *weightsFactor = forest->weightsFactor;
1163c7eeac06SToby Isaac   PetscFunctionReturn(0);
1164c7eeac06SToby Isaac }
1165c7eeac06SToby Isaac 
11669be51f97SToby Isaac /*@
11679be51f97SToby Isaac   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
11689be51f97SToby Isaac 
11699be51f97SToby Isaac   Not collective
11709be51f97SToby Isaac 
11719be51f97SToby Isaac   Input Parameter:
11729be51f97SToby Isaac . dm - the forest
11739be51f97SToby Isaac 
11749be51f97SToby Isaac   Output Parameters:
11759be51f97SToby Isaac + cStart - the first cell on this process
11769be51f97SToby Isaac - cEnd - one after the final cell on this process
11779be51f97SToby Isaac 
11781a244344SSatish Balay   Level: intermediate
11799be51f97SToby Isaac 
1180db781477SPatrick Sanan .seealso: `DMForestGetCellSF()`
11819be51f97SToby Isaac @*/
1182d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1183d71ae5a4SJacob Faibussowitsch {
1184c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1185c7eeac06SToby Isaac 
1186c7eeac06SToby Isaac   PetscFunctionBegin;
1187c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1188c7eeac06SToby Isaac   PetscValidIntPointer(cStart, 2);
1189064a246eSJacob Faibussowitsch   PetscValidIntPointer(cEnd, 3);
119048a46eb9SPierre Jolivet   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) PetscCall(forest->createcellchart(dm, &forest->cStart, &forest->cEnd));
1191c7eeac06SToby Isaac   *cStart = forest->cStart;
1192c7eeac06SToby Isaac   *cEnd   = forest->cEnd;
1193c7eeac06SToby Isaac   PetscFunctionReturn(0);
1194c7eeac06SToby Isaac }
1195c7eeac06SToby Isaac 
11969be51f97SToby Isaac /*@
11979be51f97SToby Isaac   DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes
11989be51f97SToby Isaac 
11999be51f97SToby Isaac   Not collective
12009be51f97SToby Isaac 
12019be51f97SToby Isaac   Input Parameter:
12029be51f97SToby Isaac . dm - the forest
12039be51f97SToby Isaac 
12049be51f97SToby Isaac   Output Parameter:
12059be51f97SToby Isaac . cellSF - the PetscSF
12069be51f97SToby Isaac 
12071a244344SSatish Balay   Level: intermediate
12089be51f97SToby Isaac 
1209db781477SPatrick Sanan .seealso: `DMForestGetCellChart()`
12109be51f97SToby Isaac @*/
1211d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1212d71ae5a4SJacob Faibussowitsch {
1213c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1214c7eeac06SToby Isaac 
1215c7eeac06SToby Isaac   PetscFunctionBegin;
1216c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1217c7eeac06SToby Isaac   PetscValidPointer(cellSF, 2);
121848a46eb9SPierre Jolivet   if ((!forest->cellSF) && forest->createcellsf) PetscCall(forest->createcellsf(dm, &forest->cellSF));
1219c7eeac06SToby Isaac   *cellSF = forest->cellSF;
1220c7eeac06SToby Isaac   PetscFunctionReturn(0);
1221c7eeac06SToby Isaac }
1222c7eeac06SToby Isaac 
12239be51f97SToby Isaac /*@C
12249be51f97SToby Isaac   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
12259be51f97SToby Isaac   DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination).  The
1226cd3c525cSToby Isaac   interpretation of the label values is up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP,
1227cd3c525cSToby Isaac   DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have been reserved as choices that should be accepted by all subtypes.
12289be51f97SToby Isaac 
12299be51f97SToby Isaac   Logically collective on dm
12309be51f97SToby Isaac 
12319be51f97SToby Isaac   Input Parameters:
12329be51f97SToby Isaac - dm - the forest
1233a1b0c543SToby Isaac + adaptLabel - the label in the pre-adaptation forest
12349be51f97SToby Isaac 
12359be51f97SToby Isaac   Level: intermediate
12369be51f97SToby Isaac 
1237db781477SPatrick Sanan .seealso `DMForestGetAdaptivityLabel()`
12389be51f97SToby Isaac @*/
1239d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel)
1240d71ae5a4SJacob Faibussowitsch {
1241c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1242c7eeac06SToby Isaac 
1243c7eeac06SToby Isaac   PetscFunctionBegin;
1244c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12451fd50544SStefano Zampini   if (adaptLabel) PetscValidHeaderSpecific(adaptLabel, DMLABEL_CLASSID, 2);
12469566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)adaptLabel));
12479566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&forest->adaptLabel));
12481fd50544SStefano Zampini   forest->adaptLabel = adaptLabel;
1249c7eeac06SToby Isaac   PetscFunctionReturn(0);
1250c7eeac06SToby Isaac }
1251c7eeac06SToby Isaac 
12529be51f97SToby Isaac /*@C
12539be51f97SToby Isaac   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that
12549be51f97SToby Isaac   holds the adaptation flags (refinement, coarsening, or some combination).  The interpretation of the label values is
1255cd3c525cSToby Isaac   up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have
1256cd3c525cSToby Isaac   been reserved as choices that should be accepted by all subtypes.
12579be51f97SToby Isaac 
12589be51f97SToby Isaac   Not collective
12599be51f97SToby Isaac 
12609be51f97SToby Isaac   Input Parameter:
12619be51f97SToby Isaac . dm - the forest
12629be51f97SToby Isaac 
12639be51f97SToby Isaac   Output Parameter:
12649be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest
12659be51f97SToby Isaac 
12669be51f97SToby Isaac   Level: intermediate
12679be51f97SToby Isaac 
1268db781477SPatrick Sanan .seealso `DMForestSetAdaptivityLabel()`
12699be51f97SToby Isaac @*/
1270d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel)
1271d71ae5a4SJacob Faibussowitsch {
1272c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1273c7eeac06SToby Isaac 
1274c7eeac06SToby Isaac   PetscFunctionBegin;
1275c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1276ba936b91SToby Isaac   *adaptLabel = forest->adaptLabel;
1277c7eeac06SToby Isaac   PetscFunctionReturn(0);
1278c7eeac06SToby Isaac }
1279c7eeac06SToby Isaac 
12809be51f97SToby Isaac /*@
12819be51f97SToby Isaac   DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
12829be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
12839be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
12849be51f97SToby Isaac   of 1.
12859be51f97SToby Isaac 
12869be51f97SToby Isaac   Logically collective on dm
12879be51f97SToby Isaac 
12889be51f97SToby Isaac   Input Parameters:
12899be51f97SToby Isaac + dm - the forest
12909be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
12919be51f97SToby Isaac - copyMode - how weights should reference weights
12929be51f97SToby Isaac 
12939be51f97SToby Isaac   Level: advanced
12949be51f97SToby Isaac 
1295db781477SPatrick Sanan .seealso: `DMForestGetCellWeights()`, `DMForestSetWeightCapacity()`
12969be51f97SToby Isaac @*/
1297d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1298d71ae5a4SJacob Faibussowitsch {
1299c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1300c7eeac06SToby Isaac   PetscInt   cStart, cEnd;
1301c7eeac06SToby Isaac 
1302c7eeac06SToby Isaac   PetscFunctionBegin;
1303c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
13049566063dSJacob Faibussowitsch   PetscCall(DMForestGetCellChart(dm, &cStart, &cEnd));
130563a3b9bcSJacob Faibussowitsch   PetscCheck(cEnd >= cStart, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "cell chart [%" PetscInt_FMT ",%" PetscInt_FMT ") is not valid", cStart, cEnd);
1306c7eeac06SToby Isaac   if (copyMode == PETSC_COPY_VALUES) {
130748a46eb9SPierre Jolivet     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) PetscCall(PetscMalloc1(cEnd - cStart, &forest->cellWeights));
13089566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(forest->cellWeights, weights, cEnd - cStart));
1309c7eeac06SToby Isaac     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
1310c7eeac06SToby Isaac     PetscFunctionReturn(0);
1311c7eeac06SToby Isaac   }
131248a46eb9SPierre Jolivet   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) PetscCall(PetscFree(forest->cellWeights));
1313c7eeac06SToby Isaac   forest->cellWeights         = weights;
1314c7eeac06SToby Isaac   forest->cellWeightsCopyMode = copyMode;
1315c7eeac06SToby Isaac   PetscFunctionReturn(0);
1316c7eeac06SToby Isaac }
1317c7eeac06SToby Isaac 
13189be51f97SToby Isaac /*@
13199be51f97SToby Isaac   DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
13209be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
13219be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
13229be51f97SToby Isaac   of 1.
13239be51f97SToby Isaac 
13249be51f97SToby Isaac   Not collective
13259be51f97SToby Isaac 
13269be51f97SToby Isaac   Input Parameter:
13279be51f97SToby Isaac . dm - the forest
13289be51f97SToby Isaac 
13299be51f97SToby Isaac   Output Parameter:
13309be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
13319be51f97SToby Isaac 
13329be51f97SToby Isaac   Level: advanced
13339be51f97SToby Isaac 
1334db781477SPatrick Sanan .seealso: `DMForestSetCellWeights()`, `DMForestSetWeightCapacity()`
13359be51f97SToby Isaac @*/
1336d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1337d71ae5a4SJacob Faibussowitsch {
1338c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1339c7eeac06SToby Isaac 
1340c7eeac06SToby Isaac   PetscFunctionBegin;
1341c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1342c7eeac06SToby Isaac   PetscValidPointer(weights, 2);
1343c7eeac06SToby Isaac   *weights = forest->cellWeights;
1344c7eeac06SToby Isaac   PetscFunctionReturn(0);
1345c7eeac06SToby Isaac }
1346c7eeac06SToby Isaac 
13479be51f97SToby Isaac /*@
13489be51f97SToby Isaac   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
13499be51f97SToby Isaac   a pre-adaptation forest (see DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each
13509be51f97SToby Isaac   process's cells to the process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that
13519be51f97SToby Isaac   the current process should not have any cells after repartitioning.
13529be51f97SToby Isaac 
13539be51f97SToby Isaac   Logically Collective on dm
13549be51f97SToby Isaac 
13559be51f97SToby Isaac   Input parameters:
13569be51f97SToby Isaac + dm - the forest
13579be51f97SToby Isaac - capacity - this process's capacity
13589be51f97SToby Isaac 
13599be51f97SToby Isaac   Level: advanced
13609be51f97SToby Isaac 
1361db781477SPatrick Sanan .seealso `DMForestGetWeightCapacity()`, `DMForestSetCellWeights()`, `DMForestSetCellWeightFactor()`
13629be51f97SToby Isaac @*/
1363d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1364d71ae5a4SJacob Faibussowitsch {
1365c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1366c7eeac06SToby Isaac 
1367c7eeac06SToby Isaac   PetscFunctionBegin;
1368c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
136928b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the weight capacity after setup");
137063a3b9bcSJacob Faibussowitsch   PetscCheck(capacity >= 0., PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot have negative weight capacity; %g", (double)capacity);
1371c7eeac06SToby Isaac   forest->weightCapacity = capacity;
1372c7eeac06SToby Isaac   PetscFunctionReturn(0);
1373c7eeac06SToby Isaac }
1374c7eeac06SToby Isaac 
13759be51f97SToby Isaac /*@
13769be51f97SToby Isaac   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
13779be51f97SToby Isaac   DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each process's cells to the
13789be51f97SToby Isaac   process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that the current process
13799be51f97SToby Isaac   should not have any cells after repartitioning.
13809be51f97SToby Isaac 
13819be51f97SToby Isaac   Not collective
13829be51f97SToby Isaac 
13839be51f97SToby Isaac   Input parameter:
13849be51f97SToby Isaac . dm - the forest
13859be51f97SToby Isaac 
13869be51f97SToby Isaac   Output parameter:
13879be51f97SToby Isaac . capacity - this process's capacity
13889be51f97SToby Isaac 
13899be51f97SToby Isaac   Level: advanced
13909be51f97SToby Isaac 
1391db781477SPatrick Sanan .seealso `DMForestSetWeightCapacity()`, `DMForestSetCellWeights()`, `DMForestSetCellWeightFactor()`
13929be51f97SToby Isaac @*/
1393d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1394d71ae5a4SJacob Faibussowitsch {
1395c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1396c7eeac06SToby Isaac 
1397c7eeac06SToby Isaac   PetscFunctionBegin;
1398c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1399c7eeac06SToby Isaac   PetscValidRealPointer(capacity, 2);
1400c7eeac06SToby Isaac   *capacity = forest->weightCapacity;
1401c7eeac06SToby Isaac   PetscFunctionReturn(0);
1402c7eeac06SToby Isaac }
1403c7eeac06SToby Isaac 
1404d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(DM dm, PetscOptionItems *PetscOptionsObject)
1405d71ae5a4SJacob Faibussowitsch {
140656ba9f64SToby Isaac   PetscBool                  flg, flg1, flg2, flg3, flg4;
1407dd8e54a2SToby Isaac   DMForestTopology           oldTopo;
1408c7eeac06SToby Isaac   char                       stringBuffer[256];
1409dd8e54a2SToby Isaac   PetscViewer                viewer;
1410dd8e54a2SToby Isaac   PetscViewerFormat          format;
141156ba9f64SToby Isaac   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1412c7eeac06SToby Isaac   PetscReal                  weightsFactor;
1413c7eeac06SToby Isaac   DMForestAdaptivityStrategy adaptStrategy;
1414db4d5e8cSToby Isaac 
1415db4d5e8cSToby Isaac   PetscFunctionBegin;
14169566063dSJacob Faibussowitsch   PetscCall(DMForestGetTopology(dm, &oldTopo));
1417d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "DMForest Options");
14189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-dm_forest_topology", "the topology of the forest's base mesh", "DMForestSetTopology", oldTopo, stringBuffer, sizeof(stringBuffer), &flg1));
14199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsViewer("-dm_forest_base_dm", "load the base DM from a viewer specification", "DMForestSetBaseDM", &viewer, &format, &flg2));
14209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsViewer("-dm_forest_coarse_forest", "load the coarse forest from a viewer specification", "DMForestSetCoarseForest", &viewer, &format, &flg3));
14219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsViewer("-dm_forest_fine_forest", "load the fine forest from a viewer specification", "DMForestSetFineForest", &viewer, &format, &flg4));
14221dca8a05SBarry Smith   PetscCheck((PetscInt)flg1 + (PetscInt)flg2 + (PetscInt)flg3 + (PetscInt)flg4 <= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "Specify only one of -dm_forest_{topology,base_dm,coarse_forest,fine_forest}");
142356ba9f64SToby Isaac   if (flg1) {
14249566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm, (DMForestTopology)stringBuffer));
14259566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm, NULL));
14269566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, NULL));
142756ba9f64SToby Isaac   }
142856ba9f64SToby Isaac   if (flg2) {
1429dd8e54a2SToby Isaac     DM base;
1430dd8e54a2SToby Isaac 
14319566063dSJacob Faibussowitsch     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &base));
14329566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
14339566063dSJacob Faibussowitsch     PetscCall(DMLoad(base, viewer));
14349566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
14359566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm, base));
14369566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&base));
14379566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm, NULL));
14389566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, NULL));
1439dd8e54a2SToby Isaac   }
144056ba9f64SToby Isaac   if (flg3) {
1441dd8e54a2SToby Isaac     DM coarse;
1442dd8e54a2SToby Isaac 
14439566063dSJacob Faibussowitsch     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &coarse));
14449566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
14459566063dSJacob Faibussowitsch     PetscCall(DMLoad(coarse, viewer));
14469566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
14479566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, coarse));
14489566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&coarse));
14499566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm, NULL));
14509566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm, NULL));
1451dd8e54a2SToby Isaac   }
145256ba9f64SToby Isaac   if (flg4) {
1453dd8e54a2SToby Isaac     DM fine;
1454dd8e54a2SToby Isaac 
14559566063dSJacob Faibussowitsch     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &fine));
14569566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
14579566063dSJacob Faibussowitsch     PetscCall(DMLoad(fine, viewer));
14589566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
14599566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, fine));
14609566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&fine));
14619566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm, NULL));
14629566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm, NULL));
1463dd8e54a2SToby Isaac   }
14649566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdjacencyDimension(dm, &adjDim));
14659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_dimension", "set the dimension of points that define adjacency in the forest", "DMForestSetAdjacencyDimension", adjDim, &adjDim, &flg, 0));
1466dd8e54a2SToby Isaac   if (flg) {
14679566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdjacencyDimension(dm, adjDim));
1468f885a11aSToby Isaac   } else {
14699566063dSJacob Faibussowitsch     PetscCall(DMForestGetAdjacencyCodimension(dm, &adjCodim));
14709566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_codimension", "set the codimension of points that define adjacency in the forest", "DMForestSetAdjacencyCodimension", adjCodim, &adjCodim, &flg, 1));
14711baa6e33SBarry Smith     if (flg) PetscCall(DMForestSetAdjacencyCodimension(dm, adjCodim));
1472dd8e54a2SToby Isaac   }
14739566063dSJacob Faibussowitsch   PetscCall(DMForestGetPartitionOverlap(dm, &overlap));
14749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_partition_overlap", "set the degree of partition overlap", "DMForestSetPartitionOverlap", overlap, &overlap, &flg, 0));
14751baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetPartitionOverlap(dm, overlap));
1476a6121fbdSMatthew G. Knepley #if 0
14779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_refine","equivalent to -dm_forest_set_minimum_refinement and -dm_forest_set_initial_refinement with the same value",NULL,minRefinement,&minRefinement,&flg,0));
1478a6121fbdSMatthew G. Knepley   if (flg) {
14799566063dSJacob Faibussowitsch     PetscCall(DMForestSetMinimumRefinement(dm,minRefinement));
14809566063dSJacob Faibussowitsch     PetscCall(DMForestSetInitialRefinement(dm,minRefinement));
1481a6121fbdSMatthew G. Knepley   }
14829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_refine_hierarchy","equivalent to -dm_forest_set_minimum_refinement 0 and -dm_forest_set_initial_refinement",NULL,initRefinement,&initRefinement,&flg,0));
1483a6121fbdSMatthew G. Knepley   if (flg) {
14849566063dSJacob Faibussowitsch     PetscCall(DMForestSetMinimumRefinement(dm,0));
14859566063dSJacob Faibussowitsch     PetscCall(DMForestSetInitialRefinement(dm,initRefinement));
1486a6121fbdSMatthew G. Knepley   }
1487a6121fbdSMatthew G. Knepley #endif
14889566063dSJacob Faibussowitsch   PetscCall(DMForestGetMinimumRefinement(dm, &minRefinement));
14899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_minimum_refinement", "set the minimum level of refinement in the forest", "DMForestSetMinimumRefinement", minRefinement, &minRefinement, &flg, 0));
14901baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetMinimumRefinement(dm, minRefinement));
14919566063dSJacob Faibussowitsch   PetscCall(DMForestGetInitialRefinement(dm, &initRefinement));
14929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_initial_refinement", "set the initial level of refinement in the forest", "DMForestSetInitialRefinement", initRefinement, &initRefinement, &flg, 0));
14931baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetInitialRefinement(dm, initRefinement));
14949566063dSJacob Faibussowitsch   PetscCall(DMForestGetMaximumRefinement(dm, &maxRefinement));
14959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_maximum_refinement", "set the maximum level of refinement in the forest", "DMForestSetMaximumRefinement", maxRefinement, &maxRefinement, &flg, 0));
14961baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetMaximumRefinement(dm, maxRefinement));
14979566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivityStrategy(dm, &adaptStrategy));
14989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-dm_forest_adaptivity_strategy", "the forest's adaptivity-flag resolution strategy", "DMForestSetAdaptivityStrategy", adaptStrategy, stringBuffer, sizeof(stringBuffer), &flg));
14991baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetAdaptivityStrategy(dm, (DMForestAdaptivityStrategy)stringBuffer));
15009566063dSJacob Faibussowitsch   PetscCall(DMForestGetGradeFactor(dm, &grade));
15019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_grade_factor", "grade factor between neighboring cells", "DMForestSetGradeFactor", grade, &grade, &flg, 0));
15021baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetGradeFactor(dm, grade));
15039566063dSJacob Faibussowitsch   PetscCall(DMForestGetCellWeightFactor(dm, &weightsFactor));
15049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_forest_cell_weight_factor", "multiplying weight factor for cell refinement", "DMForestSetCellWeightFactor", weightsFactor, &weightsFactor, &flg));
15051baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetCellWeightFactor(dm, weightsFactor));
1506d0609cedSBarry Smith   PetscOptionsHeadEnd();
1507db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1508db4d5e8cSToby Isaac }
1509db4d5e8cSToby Isaac 
1510d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1511d71ae5a4SJacob Faibussowitsch {
1512d8984e3bSMatthew G. Knepley   PetscFunctionBegin;
15139566063dSJacob Faibussowitsch   if (subdm) PetscCall(DMClone(dm, subdm));
15149566063dSJacob Faibussowitsch   PetscCall(DMCreateSectionSubDM(dm, numFields, fields, is, subdm));
1515d8984e3bSMatthew G. Knepley   PetscFunctionReturn(0);
1516d8984e3bSMatthew G. Knepley }
1517d8984e3bSMatthew G. Knepley 
1518d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
1519d71ae5a4SJacob Faibussowitsch {
15205421bac9SToby Isaac   DMLabel refine;
15215421bac9SToby Isaac   DM      fineDM;
15225421bac9SToby Isaac 
15235421bac9SToby Isaac   PetscFunctionBegin;
15249566063dSJacob Faibussowitsch   PetscCall(DMGetFineDM(dm, &fineDM));
15255421bac9SToby Isaac   if (fineDM) {
15269566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fineDM));
15275421bac9SToby Isaac     *dmRefined = fineDM;
15285421bac9SToby Isaac     PetscFunctionReturn(0);
15295421bac9SToby Isaac   }
15309566063dSJacob Faibussowitsch   PetscCall(DMForestTemplate(dm, comm, dmRefined));
15319566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "refine", &refine));
15325421bac9SToby Isaac   if (!refine) {
15339566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "refine", &refine));
15349566063dSJacob Faibussowitsch     PetscCall(DMLabelSetDefaultValue(refine, DM_ADAPT_REFINE));
15351baa6e33SBarry Smith   } else PetscCall(PetscObjectReference((PetscObject)refine));
15369566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityLabel(*dmRefined, refine));
15379566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&refine));
15385421bac9SToby Isaac   PetscFunctionReturn(0);
15395421bac9SToby Isaac }
15405421bac9SToby Isaac 
1541d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
1542d71ae5a4SJacob Faibussowitsch {
15435421bac9SToby Isaac   DMLabel coarsen;
15445421bac9SToby Isaac   DM      coarseDM;
15455421bac9SToby Isaac 
15465421bac9SToby Isaac   PetscFunctionBegin;
15474098eed7SToby Isaac   {
15484098eed7SToby Isaac     PetscMPIInt mpiComparison;
15494098eed7SToby Isaac     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);
15504098eed7SToby Isaac 
15519566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_compare(comm, dmcomm, &mpiComparison));
15521dca8a05SBarry Smith     PetscCheck(mpiComparison == MPI_IDENT || mpiComparison == MPI_CONGRUENT, dmcomm, PETSC_ERR_SUP, "No support for different communicators yet");
15534098eed7SToby Isaac   }
15549566063dSJacob Faibussowitsch   PetscCall(DMGetCoarseDM(dm, &coarseDM));
15555421bac9SToby Isaac   if (coarseDM) {
15569566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)coarseDM));
15575421bac9SToby Isaac     *dmCoarsened = coarseDM;
15585421bac9SToby Isaac     PetscFunctionReturn(0);
15595421bac9SToby Isaac   }
15609566063dSJacob Faibussowitsch   PetscCall(DMForestTemplate(dm, comm, dmCoarsened));
15619566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityPurpose(*dmCoarsened, DM_ADAPT_COARSEN));
15629566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "coarsen", &coarsen));
15635421bac9SToby Isaac   if (!coarsen) {
15649566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "coarsen", &coarsen));
15659566063dSJacob Faibussowitsch     PetscCall(DMLabelSetDefaultValue(coarsen, DM_ADAPT_COARSEN));
15661baa6e33SBarry Smith   } else PetscCall(PetscObjectReference((PetscObject)coarsen));
15679566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityLabel(*dmCoarsened, coarsen));
15689566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&coarsen));
15695421bac9SToby Isaac   PetscFunctionReturn(0);
15705421bac9SToby Isaac }
15715421bac9SToby Isaac 
1572d71ae5a4SJacob Faibussowitsch PetscErrorCode DMAdaptLabel_Forest(DM dm, PETSC_UNUSED Vec metric, DMLabel label, PETSC_UNUSED DMLabel rgLabel, DM *adaptedDM)
1573d71ae5a4SJacob Faibussowitsch {
157409350103SToby Isaac   PetscBool success;
157509350103SToby Isaac 
157609350103SToby Isaac   PetscFunctionBegin;
15779566063dSJacob Faibussowitsch   PetscCall(DMForestTemplate(dm, PetscObjectComm((PetscObject)dm), adaptedDM));
15789566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityLabel(*adaptedDM, label));
15799566063dSJacob Faibussowitsch   PetscCall(DMSetUp(*adaptedDM));
15809566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivitySuccess(*adaptedDM, &success));
158109350103SToby Isaac   if (!success) {
15829566063dSJacob Faibussowitsch     PetscCall(DMDestroy(adaptedDM));
158309350103SToby Isaac     *adaptedDM = NULL;
158409350103SToby Isaac   }
158509350103SToby Isaac   PetscFunctionReturn(0);
158609350103SToby Isaac }
158709350103SToby Isaac 
1588d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMInitialize_Forest(DM dm)
1589d71ae5a4SJacob Faibussowitsch {
1590d222f98bSToby Isaac   PetscFunctionBegin;
15919566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(dm->ops, sizeof(*(dm->ops))));
1592d222f98bSToby Isaac 
1593d222f98bSToby Isaac   dm->ops->clone          = DMClone_Forest;
1594d222f98bSToby Isaac   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1595d222f98bSToby Isaac   dm->ops->destroy        = DMDestroy_Forest;
1596d8984e3bSMatthew G. Knepley   dm->ops->createsubdm    = DMCreateSubDM_Forest;
15975421bac9SToby Isaac   dm->ops->refine         = DMRefine_Forest;
15985421bac9SToby Isaac   dm->ops->coarsen        = DMCoarsen_Forest;
1599d222f98bSToby Isaac   PetscFunctionReturn(0);
1600d222f98bSToby Isaac }
1601d222f98bSToby Isaac 
16029be51f97SToby Isaac /*MC
16039be51f97SToby Isaac 
1604bae1f979SBarry Smith      DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh.  Forests usually have a base DM
1605bae1f979SBarry Smith   (see DMForestGetBaseDM()), from which it is refined.  The refinement and partitioning of forests is considered
1606bae1f979SBarry Smith   immutable after DMSetUp() is called.  To adapt a mesh, one should call DMForestTemplate() to create a new mesh that
1607bae1f979SBarry Smith   will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new
1608bae1f979SBarry Smith   mesh.
1609bae1f979SBarry Smith 
1610bae1f979SBarry Smith   To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the
1611bae1f979SBarry Smith   previous mesh whose values indicate which cells should be refined (DM_ADAPT_REFINE) or coarsened (DM_ADAPT_COARSEN)
1612bae1f979SBarry Smith   and how (subtypes are free to allow additional values for things like anisotropic refinement).  The label should be
1613bae1f979SBarry Smith   given to the *new* mesh with DMForestSetAdaptivityLabel().
16149be51f97SToby Isaac 
16159be51f97SToby Isaac   Level: advanced
16169be51f97SToby Isaac 
1617db781477SPatrick Sanan .seealso: `DMType`, `DMCreate()`, `DMSetType()`, `DMForestGetBaseDM()`, `DMForestSetBaseDM()`, `DMForestTemplate()`, `DMForestSetAdaptivityLabel()`
16189be51f97SToby Isaac M*/
16199be51f97SToby Isaac 
1620d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1621d71ae5a4SJacob Faibussowitsch {
1622db4d5e8cSToby Isaac   DM_Forest *forest;
1623db4d5e8cSToby Isaac 
1624db4d5e8cSToby Isaac   PetscFunctionBegin;
1625db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
16264dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&forest));
1627db4d5e8cSToby Isaac   dm->dim                     = 0;
1628db4d5e8cSToby Isaac   dm->data                    = forest;
1629db4d5e8cSToby Isaac   forest->refct               = 1;
1630db4d5e8cSToby Isaac   forest->data                = NULL;
1631db4d5e8cSToby Isaac   forest->topology            = NULL;
1632cd3c525cSToby Isaac   forest->adapt               = NULL;
1633db4d5e8cSToby Isaac   forest->base                = NULL;
16346a87ffbfSToby Isaac   forest->adaptPurpose        = DM_ADAPT_DETERMINE;
1635db4d5e8cSToby Isaac   forest->adjDim              = PETSC_DEFAULT;
1636db4d5e8cSToby Isaac   forest->overlap             = PETSC_DEFAULT;
1637db4d5e8cSToby Isaac   forest->minRefinement       = PETSC_DEFAULT;
1638db4d5e8cSToby Isaac   forest->maxRefinement       = PETSC_DEFAULT;
163956ba9f64SToby Isaac   forest->initRefinement      = PETSC_DEFAULT;
1640c7eeac06SToby Isaac   forest->cStart              = PETSC_DETERMINE;
1641c7eeac06SToby Isaac   forest->cEnd                = PETSC_DETERMINE;
1642cd3c525cSToby Isaac   forest->cellSF              = NULL;
1643ebdf65a2SToby Isaac   forest->adaptLabel          = NULL;
1644db4d5e8cSToby Isaac   forest->gradeFactor         = 2;
1645db4d5e8cSToby Isaac   forest->cellWeights         = NULL;
1646db4d5e8cSToby Isaac   forest->cellWeightsCopyMode = PETSC_USE_POINTER;
1647db4d5e8cSToby Isaac   forest->weightsFactor       = 1.;
1648db4d5e8cSToby Isaac   forest->weightCapacity      = 1.;
16499566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityStrategy(dm, DMFORESTADAPTALL));
16509566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Forest(dm));
1651db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1652db4d5e8cSToby Isaac }
1653