xref: /petsc/src/dm/impls/forest/forest.c (revision 60225df5d8469840be2bf9c1f64795a92b19f3c2)
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   }
283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2927d4645fSToby Isaac }
3027d4645fSToby Isaac 
31d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMForestPackageInitialize(void)
32d71ae5a4SJacob Faibussowitsch {
3327d4645fSToby Isaac   PetscFunctionBegin;
343ba16761SJacob Faibussowitsch   if (DMForestPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
3527d4645fSToby Isaac   DMForestPackageInitialized = PETSC_TRUE;
36f885a11aSToby Isaac 
379566063dSJacob Faibussowitsch   PetscCall(DMForestRegisterType(DMFOREST));
389566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(DMForestPackageFinalize));
393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4027d4645fSToby Isaac }
4127d4645fSToby Isaac 
429be51f97SToby Isaac /*@C
43dce8aebaSBarry Smith   DMForestRegisterType - Registers a `DMType` as a subtype of `DMFOREST` (so that `DMIsForest()` will be correct)
449be51f97SToby Isaac 
459be51f97SToby Isaac   Not Collective
469be51f97SToby Isaac 
47*60225df5SJacob Faibussowitsch   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;
643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 
72*60225df5SJacob Faibussowitsch   Input Parameter:
739be51f97SToby Isaac . dm - the DM object
749be51f97SToby Isaac 
75*60225df5SJacob Faibussowitsch   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;
923ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
9327d4645fSToby Isaac     }
9427d4645fSToby Isaac     link = link->next;
9527d4645fSToby Isaac   }
9627d4645fSToby Isaac   *isForest = PETSC_FALSE;
973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9827d4645fSToby Isaac }
9927d4645fSToby Isaac 
1009be51f97SToby Isaac /*@
10120f4b53cSBarry Smith   DMForestTemplate - Create a new `DM` that will be adapted from a source `DM`.  The new `DM` reproduces the configuration
10220f4b53cSBarry Smith   of the source, but is not yet setup, so that the user can then define only the ways that the new `DM` should differ
10320f4b53cSBarry Smith   (by, e.g., refinement or repartitioning).  The source `DM` is also set as the adaptivity source `DM` of the new `DM` (see
10420f4b53cSBarry Smith   `DMForestSetAdaptivityForest()`).
1059be51f97SToby Isaac 
10620f4b53cSBarry Smith   Collective
1079be51f97SToby Isaac 
1089be51f97SToby Isaac   Input Parameters:
10920f4b53cSBarry Smith + dm   - the source `DM` object
11020f4b53cSBarry Smith - 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:
11320f4b53cSBarry Smith . tdm - the new `DM` object
1149be51f97SToby Isaac 
1159be51f97SToby Isaac   Level: intermediate
1169be51f97SToby Isaac 
11720f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `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));
1683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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));
1843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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;
1923ba16761SJacob Faibussowitsch   if (--forest->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
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));
2033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
204db4d5e8cSToby Isaac }
205db4d5e8cSToby Isaac 
2069be51f97SToby Isaac /*@C
20720f4b53cSBarry Smith   DMForestSetTopology - Set the topology of a `DMFOREST` during the pre-setup phase.  The topology is a string (e.g.
20820f4b53cSBarry Smith   "cube", "shell") and can be interpreted by subtypes of `DMFOREST`) to construct the base DM of a forest during
20920f4b53cSBarry Smith   `DMSetUp()`.
2109be51f97SToby Isaac 
21120f4b53cSBarry Smith   Logically collectiv
2129be51f97SToby Isaac 
213*60225df5SJacob Faibussowitsch   Input Parameters:
2149be51f97SToby Isaac + dm       - the forest
2159be51f97SToby Isaac - topology - the topology of the forest
2169be51f97SToby Isaac 
2179be51f97SToby Isaac   Level: intermediate
2189be51f97SToby Isaac 
21920f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `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));
2303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
231db4d5e8cSToby Isaac }
232db4d5e8cSToby Isaac 
2339be51f97SToby Isaac /*@C
23420f4b53cSBarry Smith   DMForestGetTopology - Get a string describing the topology of a `DMFOREST`.
2359be51f97SToby Isaac 
23620f4b53cSBarry Smith   Not Collective
2379be51f97SToby Isaac 
238*60225df5SJacob Faibussowitsch   Input Parameter:
2399be51f97SToby Isaac . dm - the forest
2409be51f97SToby Isaac 
241*60225df5SJacob Faibussowitsch   Output Parameter:
2429be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell')
2439be51f97SToby Isaac 
2449be51f97SToby Isaac   Level: intermediate
2459be51f97SToby Isaac 
24620f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `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;
2563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
257dd8e54a2SToby Isaac }
258dd8e54a2SToby Isaac 
2599be51f97SToby Isaac /*@
26020f4b53cSBarry Smith   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 
26420f4b53cSBarry Smith   Logically Collective
2659be51f97SToby Isaac 
2669be51f97SToby Isaac   Input Parameters:
2679be51f97SToby Isaac + dm   - the forest
26820f4b53cSBarry Smith - base - the base `DM` of the forest
269765b024eSBarry Smith 
2709be51f97SToby Isaac   Level: intermediate
2719be51f97SToby Isaac 
27220f4b53cSBarry Smith   Note:
27320f4b53cSBarry Smith   Currently the base `DM` must be a `DMPLEX`
27420f4b53cSBarry Smith 
27520f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `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));
2993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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 
30720f4b53cSBarry Smith   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 
32020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `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;
3303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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;
3413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
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;
3523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
353cf38a08cSToby Isaac }
354cf38a08cSToby Isaac 
3559be51f97SToby Isaac /*@
3569be51f97SToby Isaac   DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be
35720f4b53cSBarry Smith   adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) in `DMSetUp()`.  Usually not needed
35820f4b53cSBarry Smith   by users directly: `DMForestTemplate()` constructs a new forest to be adapted from an old forest and calls this
3599be51f97SToby Isaac   routine.
3609be51f97SToby Isaac 
36120f4b53cSBarry Smith   Logically Collective
3629be51f97SToby Isaac 
363d8d19677SJose E. Roman   Input Parameters:
3649be51f97SToby Isaac + dm    - the new forest, which will be constructed from adapt
3659be51f97SToby Isaac - adapt - the old forest
3669be51f97SToby Isaac 
3679be51f97SToby Isaac   Level: intermediate
3689be51f97SToby Isaac 
36920f4b53cSBarry Smith   Note:
37020f4b53cSBarry Smith   This can be called after setup with `adapt` = `NULL`, which will clear all internal data related to the
37120f4b53cSBarry Smith   adaptivity forest from `dm`.  This way, repeatedly adapting does not leave stale `DM` objects in memory.
37220f4b53cSBarry Smith 
37320f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdaptivityForest()`, `DMForestSetAdaptivityPurpose()`
3749be51f97SToby Isaac @*/
375d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityForest(DM dm, DM adapt)
376d71ae5a4SJacob Faibussowitsch {
377dffe73a3SToby Isaac   DM_Forest *forest, *adaptForest, *oldAdaptForest;
378dffe73a3SToby Isaac   DM         oldAdapt;
379456cc5b7SMatthew G. Knepley   PetscBool  isForest;
380dd8e54a2SToby Isaac 
381dd8e54a2SToby Isaac   PetscFunctionBegin;
382dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3831fd50544SStefano Zampini   if (adapt) PetscValidHeaderSpecific(adapt, DM_CLASSID, 2);
3849566063dSJacob Faibussowitsch   PetscCall(DMIsForest(dm, &isForest));
3853ba16761SJacob Faibussowitsch   if (!isForest) PetscFunctionReturn(PETSC_SUCCESS);
3861dca8a05SBarry Smith   PetscCheck(adapt == NULL || !dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adaptation forest after setup");
387ba936b91SToby Isaac   forest = (DM_Forest *)dm->data;
3889566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivityForest(dm, &oldAdapt));
389193eb951SToby Isaac   adaptForest    = (DM_Forest *)(adapt ? adapt->data : NULL);
390193eb951SToby Isaac   oldAdaptForest = (DM_Forest *)(oldAdapt ? oldAdapt->data : NULL);
391dffe73a3SToby Isaac   if (adaptForest != oldAdaptForest) {
3929566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&forest->preCoarseToFine));
3939566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&forest->coarseToPreFine));
3949566063dSJacob Faibussowitsch     if (forest->clearadaptivityforest) PetscCall((*forest->clearadaptivityforest)(dm));
395dffe73a3SToby Isaac   }
39626d9498aSToby Isaac   switch (forest->adaptPurpose) {
397cd3c525cSToby Isaac   case DM_ADAPT_DETERMINE:
3989566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)adapt));
3999566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&(forest->adapt)));
400ba936b91SToby Isaac     forest->adapt = adapt;
40126d9498aSToby Isaac     break;
402d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_REFINE:
403d71ae5a4SJacob Faibussowitsch     PetscCall(DMSetCoarseDM(dm, adapt));
404d71ae5a4SJacob Faibussowitsch     break;
405a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
406d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_COARSEN_LAST:
407d71ae5a4SJacob Faibussowitsch     PetscCall(DMSetFineDM(dm, adapt));
408d71ae5a4SJacob Faibussowitsch     break;
409d71ae5a4SJacob Faibussowitsch   default:
410d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid adaptivity purpose");
41126d9498aSToby Isaac   }
4123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
413dd8e54a2SToby Isaac }
414dd8e54a2SToby Isaac 
4159be51f97SToby Isaac /*@
4169be51f97SToby Isaac   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.
4179be51f97SToby Isaac 
41820f4b53cSBarry Smith   Not Collective
4199be51f97SToby Isaac 
4209be51f97SToby Isaac   Input Parameter:
4219be51f97SToby Isaac . dm - the forest
4229be51f97SToby Isaac 
4239be51f97SToby Isaac   Output Parameter:
42420f4b53cSBarry Smith . adapt - the forest from which `dm` is/was adapted
4259be51f97SToby Isaac 
4269be51f97SToby Isaac   Level: intermediate
4279be51f97SToby Isaac 
42820f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdaptivityForest()`, `DMForestSetAdaptivityPurpose()`
4299be51f97SToby Isaac @*/
430d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
431d71ae5a4SJacob Faibussowitsch {
432ba936b91SToby Isaac   DM_Forest *forest;
433dd8e54a2SToby Isaac 
434dd8e54a2SToby Isaac   PetscFunctionBegin;
435dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
436ba936b91SToby Isaac   forest = (DM_Forest *)dm->data;
43726d9498aSToby Isaac   switch (forest->adaptPurpose) {
438d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_DETERMINE:
439d71ae5a4SJacob Faibussowitsch     *adapt = forest->adapt;
440d71ae5a4SJacob Faibussowitsch     break;
441d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_REFINE:
442d71ae5a4SJacob Faibussowitsch     PetscCall(DMGetCoarseDM(dm, adapt));
443d71ae5a4SJacob Faibussowitsch     break;
444a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
445d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_COARSEN_LAST:
446d71ae5a4SJacob Faibussowitsch     PetscCall(DMGetFineDM(dm, adapt));
447d71ae5a4SJacob Faibussowitsch     break;
448d71ae5a4SJacob Faibussowitsch   default:
449d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid adaptivity purpose");
45026d9498aSToby Isaac   }
4513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45226d9498aSToby Isaac }
45326d9498aSToby Isaac 
4549be51f97SToby Isaac /*@
45520f4b53cSBarry Smith   DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current `DM` is being adapted from its
45620f4b53cSBarry Smith   source (set with `DMForestSetAdaptivityForest()`) for the purpose of refinement (`DM_ADAPT_REFINE`), coarsening
45720f4b53cSBarry Smith   (`DM_ADAPT_COARSEN`), or undefined (`DM_ADAPT_DETERMINE`).  This only matters for the purposes of reference counting:
45820f4b53cSBarry Smith   during `DMDestroy()`, cyclic references can be found between `DM`s only if the cyclic reference is due to a fine/coarse
45920f4b53cSBarry Smith   relationship (see `DMSetFineDM()`/`DMSetCoarseDM()`).  If the purpose is not refinement or coarsening, and the user does
46020f4b53cSBarry Smith   not maintain a reference to the post-adaptation forest (i.e., the one created by `DMForestTemplate()`), then this can
46120f4b53cSBarry Smith   cause a memory leak.  This method is used by subtypes of `DMFOREST` when automatically constructing mesh hierarchies.
4629be51f97SToby Isaac 
46320f4b53cSBarry Smith   Logically Collective
4649be51f97SToby Isaac 
4659be51f97SToby Isaac   Input Parameters:
4669be51f97SToby Isaac + dm      - the forest
467bf2d5fbbSStefano Zampini - purpose - the adaptivity purpose
4689be51f97SToby Isaac 
4699be51f97SToby Isaac   Level: advanced
4709be51f97SToby Isaac 
47120f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestTemplate()`, `DMForestSetAdaptivityForest()`, `DMForestGetAdaptivityForest()`, `DMAdaptFlag`
4729be51f97SToby Isaac @*/
473d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose)
474d71ae5a4SJacob Faibussowitsch {
47526d9498aSToby Isaac   DM_Forest *forest;
47626d9498aSToby Isaac 
47726d9498aSToby Isaac   PetscFunctionBegin;
47826d9498aSToby Isaac   forest = (DM_Forest *)dm->data;
47928b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adaptation forest after setup");
48026d9498aSToby Isaac   if (purpose != forest->adaptPurpose) {
48126d9498aSToby Isaac     DM adapt;
48226d9498aSToby Isaac 
4839566063dSJacob Faibussowitsch     PetscCall(DMForestGetAdaptivityForest(dm, &adapt));
4849566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)adapt));
4859566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, NULL));
486f885a11aSToby Isaac 
48726d9498aSToby Isaac     forest->adaptPurpose = purpose;
488f885a11aSToby Isaac 
4899566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, adapt));
4909566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&adapt));
49126d9498aSToby Isaac   }
4923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
493dd8e54a2SToby Isaac }
494dd8e54a2SToby Isaac 
49556c0450aSToby Isaac /*@
49620f4b53cSBarry Smith   DMForestGetAdaptivityPurpose - Get whether the current `DM` is being adapted from its source (set with
49720f4b53cSBarry Smith   `DMForestSetAdaptivityForest()`) for the purpose of refinement (`DM_ADAPT_REFINE`), coarsening (`DM_ADAPT_COARSEN`),
49820f4b53cSBarry Smith   coarsening only the last level (`DM_ADAPT_COARSEN_LAST`) or undefined (`DM_ADAPT_DETERMINE`).
49920f4b53cSBarry Smith   This only matters for the purposes of reference counting: during `DMDestroy()`, cyclic
500*60225df5SJacob Faibussowitsch 
501*60225df5SJacob Faibussowitsch   References Can Be Found Between `Dm`S Only If The Cyclic Reference Is Due To A Fine/Coarse Relationship (See
50220f4b53cSBarry Smith   `DMSetFineDM()`/`DMSetCoarseDM()`).  If the purpose is not refinement or coarsening, and the user does not maintain a
50320f4b53cSBarry Smith   reference to the post-adaptation forest (i.e., the one created by `DMForestTemplate()`), then this can cause a memory
50420f4b53cSBarry Smith   leak.  This method is used by subtypes of `DMFOREST` when automatically constructing mesh hierarchies.
50556c0450aSToby Isaac 
50620f4b53cSBarry Smith   Not Collective
50756c0450aSToby Isaac 
50856c0450aSToby Isaac   Input Parameter:
50956c0450aSToby Isaac . dm - the forest
51056c0450aSToby Isaac 
51156c0450aSToby Isaac   Output Parameter:
512bf2d5fbbSStefano Zampini . purpose - the adaptivity purpose
51356c0450aSToby Isaac 
51456c0450aSToby Isaac   Level: advanced
51556c0450aSToby Isaac 
51620f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestTemplate()`, `DMForestSetAdaptivityForest()`, `DMForestGetAdaptivityForest()`, `DMAdaptFlag`
51756c0450aSToby Isaac @*/
518d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose)
519d71ae5a4SJacob Faibussowitsch {
52056c0450aSToby Isaac   DM_Forest *forest;
52156c0450aSToby Isaac 
52256c0450aSToby Isaac   PetscFunctionBegin;
52356c0450aSToby Isaac   forest   = (DM_Forest *)dm->data;
52456c0450aSToby Isaac   *purpose = forest->adaptPurpose;
5253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
52656c0450aSToby Isaac }
52756c0450aSToby Isaac 
5289be51f97SToby Isaac /*@
5299be51f97SToby Isaac   DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
5309be51f97SToby Isaac   cell adjacency (for the purposes of partitioning and overlap).
5319be51f97SToby Isaac 
53220f4b53cSBarry Smith   Logically Collective
5339be51f97SToby Isaac 
5349be51f97SToby Isaac   Input Parameters:
5359be51f97SToby Isaac + dm     - the forest
5369be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency)
5379be51f97SToby Isaac 
5389be51f97SToby Isaac   Level: intermediate
5399be51f97SToby Isaac 
54020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`, `DMForestSetPartitionOverlap()`
5419be51f97SToby Isaac @*/
542d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
543d71ae5a4SJacob Faibussowitsch {
544dd8e54a2SToby Isaac   PetscInt   dim;
545dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
546dd8e54a2SToby Isaac 
547dd8e54a2SToby Isaac   PetscFunctionBegin;
548dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
54928b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adjacency dimension after setup");
55063a3b9bcSJacob Faibussowitsch   PetscCheck(adjDim >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "adjacency dim cannot be < 0: %" PetscInt_FMT, adjDim);
5519566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
55263a3b9bcSJacob Faibussowitsch   PetscCheck(adjDim <= dim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "adjacency dim cannot be > %" PetscInt_FMT ": %" PetscInt_FMT, dim, adjDim);
553dd8e54a2SToby Isaac   forest->adjDim = adjDim;
5543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
555dd8e54a2SToby Isaac }
556dd8e54a2SToby Isaac 
5579be51f97SToby Isaac /*@
55820f4b53cSBarry Smith   DMForestSetAdjacencyCodimension - Like `DMForestSetAdjacencyDimension()`, but specified as a co-dimension (so that,
5599be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
5609be51f97SToby Isaac 
56120f4b53cSBarry Smith   Logically Collective
5629be51f97SToby Isaac 
5639be51f97SToby Isaac   Input Parameters:
5649be51f97SToby Isaac + dm       - the forest
56520f4b53cSBarry Smith - adjCodim - default is the dimension of the forest (see `DMGetDimension()`), since this is the codimension of vertices
5669be51f97SToby Isaac 
5679be51f97SToby Isaac   Level: intermediate
5689be51f97SToby Isaac 
56920f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdjacencyCodimension()`, `DMForestSetAdjacencyDimension()`
5709be51f97SToby Isaac @*/
571d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
572d71ae5a4SJacob Faibussowitsch {
573dd8e54a2SToby Isaac   PetscInt dim;
574dd8e54a2SToby Isaac 
575dd8e54a2SToby Isaac   PetscFunctionBegin;
576dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5779566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
5789566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdjacencyDimension(dm, dim - adjCodim));
5793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
580dd8e54a2SToby Isaac }
581dd8e54a2SToby Isaac 
5829be51f97SToby Isaac /*@
5839be51f97SToby Isaac   DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
5849be51f97SToby Isaac   purposes of partitioning and overlap).
5859be51f97SToby Isaac 
58620f4b53cSBarry Smith   Not Collective
5879be51f97SToby Isaac 
5889be51f97SToby Isaac   Input Parameter:
5899be51f97SToby Isaac . dm - the forest
5909be51f97SToby Isaac 
5919be51f97SToby Isaac   Output Parameter:
5929be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency)
5939be51f97SToby Isaac 
5949be51f97SToby Isaac   Level: intermediate
5959be51f97SToby Isaac 
59620f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdjacencyDimension()`, `DMForestGetAdjacencyCodimension()`, `DMForestSetPartitionOverlap()`
5979be51f97SToby Isaac @*/
598d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
599d71ae5a4SJacob Faibussowitsch {
600dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
601dd8e54a2SToby Isaac 
602dd8e54a2SToby Isaac   PetscFunctionBegin;
603dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
604dd8e54a2SToby Isaac   PetscValidIntPointer(adjDim, 2);
605dd8e54a2SToby Isaac   *adjDim = forest->adjDim;
6063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
607dd8e54a2SToby Isaac }
608dd8e54a2SToby Isaac 
6099be51f97SToby Isaac /*@
61020f4b53cSBarry Smith   DMForestGetAdjacencyCodimension - Like `DMForestGetAdjacencyDimension()`, but specified as a co-dimension (so that,
6119be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
6129be51f97SToby Isaac 
61320f4b53cSBarry Smith   Not Collective
6149be51f97SToby Isaac 
6159be51f97SToby Isaac   Input Parameter:
6169be51f97SToby Isaac . dm - the forest
6179be51f97SToby Isaac 
6189be51f97SToby Isaac   Output Parameter:
61920f4b53cSBarry Smith . adjCodim - default isthe dimension of the forest (see `DMGetDimension()`), since this is the codimension of vertices
6209be51f97SToby Isaac 
6219be51f97SToby Isaac   Level: intermediate
6229be51f97SToby Isaac 
62320f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdjacencyCodimension()`, `DMForestGetAdjacencyDimension()`
6249be51f97SToby Isaac @*/
625d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
626d71ae5a4SJacob Faibussowitsch {
627dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
628dd8e54a2SToby Isaac   PetscInt   dim;
629dd8e54a2SToby Isaac 
630dd8e54a2SToby Isaac   PetscFunctionBegin;
631dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
632dd8e54a2SToby Isaac   PetscValidIntPointer(adjCodim, 2);
6339566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
634dd8e54a2SToby Isaac   *adjCodim = dim - forest->adjDim;
6353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
636dd8e54a2SToby Isaac }
637dd8e54a2SToby Isaac 
6389be51f97SToby Isaac /*@
6399be51f97SToby Isaac   DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel
6409be51f97SToby Isaac   partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding
6419be51f97SToby Isaac   adjacent cells
6429be51f97SToby Isaac 
64320f4b53cSBarry Smith   Logically Collective
6449be51f97SToby Isaac 
6459be51f97SToby Isaac   Input Parameters:
6469be51f97SToby Isaac + dm      - the forest
6479be51f97SToby Isaac - overlap - default 0
6489be51f97SToby Isaac 
6499be51f97SToby Isaac   Level: intermediate
6509be51f97SToby Isaac 
65120f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetPartitionOverlap()`, `DMForestSetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`
6529be51f97SToby Isaac @*/
653d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap)
654d71ae5a4SJacob Faibussowitsch {
655dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
656dd8e54a2SToby Isaac 
657dd8e54a2SToby Isaac   PetscFunctionBegin;
658dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
65928b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the overlap after setup");
66063a3b9bcSJacob Faibussowitsch   PetscCheck(overlap >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "overlap cannot be < 0: %" PetscInt_FMT, overlap);
661dd8e54a2SToby Isaac   forest->overlap = overlap;
6623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
663dd8e54a2SToby Isaac }
664dd8e54a2SToby Isaac 
6659be51f97SToby Isaac /*@
6669be51f97SToby Isaac   DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values
6679be51f97SToby Isaac   > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells
6689be51f97SToby Isaac 
66920f4b53cSBarry Smith   Not Collective
6709be51f97SToby Isaac 
6719be51f97SToby Isaac   Input Parameter:
6729be51f97SToby Isaac . dm - the forest
6739be51f97SToby Isaac 
6749be51f97SToby Isaac   Output Parameter:
6759be51f97SToby Isaac . overlap - default 0
6769be51f97SToby Isaac 
6779be51f97SToby Isaac   Level: intermediate
6789be51f97SToby Isaac 
67920f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetPartitionOverlap()`, `DMForestSetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`
6809be51f97SToby Isaac @*/
681d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap)
682d71ae5a4SJacob Faibussowitsch {
683dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
684dd8e54a2SToby Isaac 
685dd8e54a2SToby Isaac   PetscFunctionBegin;
686dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
687dd8e54a2SToby Isaac   PetscValidIntPointer(overlap, 2);
688dd8e54a2SToby Isaac   *overlap = forest->overlap;
6893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
690dd8e54a2SToby Isaac }
691dd8e54a2SToby Isaac 
6929be51f97SToby Isaac /*@
6939be51f97SToby Isaac   DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base
69420f4b53cSBarry Smith   `DM`, see `DMForestGetBaseDM()`) allowed in the forest.  If the forest is being created by coarsening a previous forest
69520f4b53cSBarry Smith   (see `DMForestGetAdaptivityForest()`) this limits the amount of coarsening.
6969be51f97SToby Isaac 
69720f4b53cSBarry Smith   Logically Collective
6989be51f97SToby Isaac 
6999be51f97SToby Isaac   Input Parameters:
7009be51f97SToby Isaac + dm            - the forest
70120f4b53cSBarry Smith - minRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`)
7029be51f97SToby Isaac 
7039be51f97SToby Isaac   Level: intermediate
7049be51f97SToby Isaac 
70520f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestSetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()`
7069be51f97SToby Isaac @*/
707d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement)
708d71ae5a4SJacob Faibussowitsch {
709dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
710dd8e54a2SToby Isaac 
711dd8e54a2SToby Isaac   PetscFunctionBegin;
712dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
71328b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the minimum refinement after setup");
714dd8e54a2SToby Isaac   forest->minRefinement = minRefinement;
7153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
716dd8e54a2SToby Isaac }
717dd8e54a2SToby Isaac 
7189be51f97SToby Isaac /*@
71920f4b53cSBarry Smith   DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base `DM`, see
72020f4b53cSBarry Smith   `DMForestGetBaseDM()`) allowed in the forest.  If the forest is being created by coarsening a previous forest (see
72120f4b53cSBarry Smith   `DMForestGetAdaptivityForest()`), this limits the amount of coarsening.
7229be51f97SToby Isaac 
72320f4b53cSBarry Smith   Not Collective
7249be51f97SToby Isaac 
7259be51f97SToby Isaac   Input Parameter:
7269be51f97SToby Isaac . dm - the forest
7279be51f97SToby Isaac 
7289be51f97SToby Isaac   Output Parameter:
72920f4b53cSBarry Smith . minRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`)
7309be51f97SToby Isaac 
7319be51f97SToby Isaac   Level: intermediate
7329be51f97SToby Isaac 
73320f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetMinimumRefinement()`, `DMForestGetMaximumRefinement()`, `DMForestGetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()`
7349be51f97SToby Isaac @*/
735d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement)
736d71ae5a4SJacob Faibussowitsch {
737dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
738dd8e54a2SToby Isaac 
739dd8e54a2SToby Isaac   PetscFunctionBegin;
740dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
741dd8e54a2SToby Isaac   PetscValidIntPointer(minRefinement, 2);
742dd8e54a2SToby Isaac   *minRefinement = forest->minRefinement;
7433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
744dd8e54a2SToby Isaac }
745dd8e54a2SToby Isaac 
7469be51f97SToby Isaac /*@
7479be51f97SToby Isaac   DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base
74820f4b53cSBarry Smith   `DM`, see `DMForestGetBaseDM()`) allowed in the forest.
7499be51f97SToby Isaac 
75020f4b53cSBarry Smith   Logically Collective
7519be51f97SToby Isaac 
7529be51f97SToby Isaac   Input Parameters:
7539be51f97SToby Isaac + dm             - the forest
754*60225df5SJacob Faibussowitsch - initRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`)
7559be51f97SToby Isaac 
7569be51f97SToby Isaac   Level: intermediate
7579be51f97SToby Isaac 
75820f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestGetBaseDM()`
7599be51f97SToby Isaac @*/
760d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement)
761d71ae5a4SJacob Faibussowitsch {
76256ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
76356ba9f64SToby Isaac 
76456ba9f64SToby Isaac   PetscFunctionBegin;
76556ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
76628b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the initial refinement after setup");
76756ba9f64SToby Isaac   forest->initRefinement = initRefinement;
7683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76956ba9f64SToby Isaac }
77056ba9f64SToby Isaac 
7719be51f97SToby Isaac /*@
77220f4b53cSBarry Smith   DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base `DM`, see
77320f4b53cSBarry Smith   `DMForestGetBaseDM()`) allowed in the forest.
7749be51f97SToby Isaac 
77520f4b53cSBarry Smith   Not Collective
7769be51f97SToby Isaac 
7779be51f97SToby Isaac   Input Parameter:
7789be51f97SToby Isaac . dm - the forest
7799be51f97SToby Isaac 
78001d2d390SJose E. Roman   Output Parameter:
78120f4b53cSBarry Smith . initRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`)
7829be51f97SToby Isaac 
7839be51f97SToby Isaac   Level: intermediate
7849be51f97SToby Isaac 
78520f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestGetBaseDM()`
7869be51f97SToby Isaac @*/
787d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement)
788d71ae5a4SJacob Faibussowitsch {
78956ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
79056ba9f64SToby Isaac 
79156ba9f64SToby Isaac   PetscFunctionBegin;
79256ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
79356ba9f64SToby Isaac   PetscValidIntPointer(initRefinement, 2);
79456ba9f64SToby Isaac   *initRefinement = forest->initRefinement;
7953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
79656ba9f64SToby Isaac }
79756ba9f64SToby Isaac 
7989be51f97SToby Isaac /*@
7999be51f97SToby Isaac   DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base
80020f4b53cSBarry Smith   `DM`, see `DMForestGetBaseDM()`) allowed in the forest.  If the forest is being created by refining a previous forest
80120f4b53cSBarry Smith   (see `DMForestGetAdaptivityForest()`), this limits the amount of refinement.
8029be51f97SToby Isaac 
80320f4b53cSBarry Smith   Logically Collective
8049be51f97SToby Isaac 
8059be51f97SToby Isaac   Input Parameters:
8069be51f97SToby Isaac + dm            - the forest
80720f4b53cSBarry Smith - maxRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`)
8089be51f97SToby Isaac 
8099be51f97SToby Isaac   Level: intermediate
8109be51f97SToby Isaac 
81120f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestSetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityDM()`
8129be51f97SToby Isaac @*/
813d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement)
814d71ae5a4SJacob Faibussowitsch {
815dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
816dd8e54a2SToby Isaac 
817dd8e54a2SToby Isaac   PetscFunctionBegin;
818dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
81928b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the maximum refinement after setup");
820c7eeac06SToby Isaac   forest->maxRefinement = maxRefinement;
8213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
822dd8e54a2SToby Isaac }
823dd8e54a2SToby Isaac 
8249be51f97SToby Isaac /*@
82520f4b53cSBarry Smith   DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base `DM`, see
82620f4b53cSBarry Smith   `DMForestGetBaseDM()`) allowed in the forest.  If the forest is being created by refining a previous forest (see
827f826b5fcSPierre Jolivet   `DMForestGetAdaptivityForest()`), this limits the amount of refinement.
8289be51f97SToby Isaac 
82920f4b53cSBarry Smith   Not Collective
8309be51f97SToby Isaac 
8319be51f97SToby Isaac   Input Parameter:
8329be51f97SToby Isaac . dm - the forest
8339be51f97SToby Isaac 
8349be51f97SToby Isaac   Output Parameter:
83520f4b53cSBarry Smith . maxRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`)
8369be51f97SToby Isaac 
8379be51f97SToby Isaac   Level: intermediate
8389be51f97SToby Isaac 
83920f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetMaximumRefinement()`, `DMForestGetMinimumRefinement()`, `DMForestGetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()`
8409be51f97SToby Isaac @*/
841d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement)
842d71ae5a4SJacob Faibussowitsch {
843dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
844dd8e54a2SToby Isaac 
845dd8e54a2SToby Isaac   PetscFunctionBegin;
846dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
847c7eeac06SToby Isaac   PetscValidIntPointer(maxRefinement, 2);
848c7eeac06SToby Isaac   *maxRefinement = forest->maxRefinement;
8493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
850dd8e54a2SToby Isaac }
851c7eeac06SToby Isaac 
8529be51f97SToby Isaac /*@C
8539be51f97SToby Isaac   DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes.
8549be51f97SToby Isaac 
85520f4b53cSBarry Smith   Logically Collective
8569be51f97SToby Isaac 
8579be51f97SToby Isaac   Input Parameters:
8589be51f97SToby Isaac + dm            - the forest
85920f4b53cSBarry Smith - adaptStrategy - default `DMFORESTADAPTALL`
8609be51f97SToby Isaac 
8619be51f97SToby Isaac   Level: advanced
8629be51f97SToby Isaac 
86320f4b53cSBarry Smith   Notes:
86420f4b53cSBarry Smith   Subtypes of `DMFOREST` may define their own strategies.  Two default strategies are `DMFORESTADAPTALL`, which indicates that all processes must agree
86520f4b53cSBarry Smith   for a refinement/coarsening flag to be valid, and `DMFORESTADAPTANY`, which indicates that only one process needs to
86620f4b53cSBarry Smith   specify refinement/coarsening.
86720f4b53cSBarry Smith 
86820f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdaptivityStrategy()`, `DMFORESTADAPTALL`, `DMFORESTADAPTANY`
8699be51f97SToby Isaac @*/
870d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy)
871d71ae5a4SJacob Faibussowitsch {
872c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
873c7eeac06SToby Isaac 
874c7eeac06SToby Isaac   PetscFunctionBegin;
875c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8769566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest->adaptStrategy));
8779566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy((const char *)adaptStrategy, (char **)&forest->adaptStrategy));
8783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
879c7eeac06SToby Isaac }
880c7eeac06SToby Isaac 
8819be51f97SToby Isaac /*@C
882*60225df5SJacob Faibussowitsch   DMForestGetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes.
8839be51f97SToby Isaac 
88420f4b53cSBarry Smith   Not Collective
8859be51f97SToby Isaac 
8869be51f97SToby Isaac   Input Parameter:
8879be51f97SToby Isaac . dm - the forest
8889be51f97SToby Isaac 
8899be51f97SToby Isaac   Output Parameter:
89020f4b53cSBarry Smith . adaptStrategy - the adaptivity strategy (default `DMFORESTADAPTALL`)
8919be51f97SToby Isaac 
8929be51f97SToby Isaac   Level: advanced
8939be51f97SToby Isaac 
89420f4b53cSBarry Smith   Note:
89520f4b53cSBarry Smith   Subtypes
89620f4b53cSBarry Smith   of `DMFOREST` may define their own strategies.  Two default strategies are `DMFORESTADAPTALL`, which indicates that all
89720f4b53cSBarry Smith   processes must agree for a refinement/coarsening flag to be valid, and `DMFORESTADAPTANY`, which indicates that only
89820f4b53cSBarry Smith   one process needs to specify refinement/coarsening.
89920f4b53cSBarry Smith 
90020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMFORESTADAPTALL`, `DMFORESTADAPTANY`, `DMForestSetAdaptivityStrategy()`
9019be51f97SToby Isaac @*/
902d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
903d71ae5a4SJacob Faibussowitsch {
904c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
905c7eeac06SToby Isaac 
906c7eeac06SToby Isaac   PetscFunctionBegin;
907c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
908c7eeac06SToby Isaac   PetscValidPointer(adaptStrategy, 2);
909c7eeac06SToby Isaac   *adaptStrategy = forest->adaptStrategy;
9103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
911c7eeac06SToby Isaac }
912c7eeac06SToby Isaac 
9132a133e43SToby Isaac /*@
9142a133e43SToby Isaac   DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning,
91520f4b53cSBarry Smith   etc.) was successful.
9162a133e43SToby Isaac 
91720f4b53cSBarry Smith   Collective
9182a133e43SToby Isaac 
9192a133e43SToby Isaac   Input Parameter:
9202a133e43SToby Isaac . dm - the post-adaptation forest
9212a133e43SToby Isaac 
9222a133e43SToby Isaac   Output Parameter:
92320f4b53cSBarry Smith . success - `PETSC_TRUE` if the post-adaptation forest is different from the pre-adaptation forest.
9242a133e43SToby Isaac 
9252a133e43SToby Isaac   Level: intermediate
9262a133e43SToby Isaac 
92720f4b53cSBarry Smith   Notes:
92820f4b53cSBarry Smith   `PETSC_FALSE` indicates that the post-adaptation forest is the same as the pre-adpatation
92920f4b53cSBarry Smith   forest.  A requested adaptation may have been unsuccessful if, for example, the requested refinement would have
93020f4b53cSBarry Smith   exceeded the maximum refinement level.
93120f4b53cSBarry Smith 
93220f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`
9332a133e43SToby Isaac @*/
934d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success)
935d71ae5a4SJacob Faibussowitsch {
9362a133e43SToby Isaac   DM_Forest *forest;
9372a133e43SToby Isaac 
9382a133e43SToby Isaac   PetscFunctionBegin;
9392a133e43SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
94028b400f6SJacob Faibussowitsch   PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() has not been called yet.");
9412a133e43SToby Isaac   forest = (DM_Forest *)dm->data;
9429566063dSJacob Faibussowitsch   PetscCall((forest->getadaptivitysuccess)(dm, success));
9433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9442a133e43SToby Isaac }
9452a133e43SToby Isaac 
946bf9b5d84SToby Isaac /*@
94720f4b53cSBarry Smith   DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer `PetscSF`s should be computed
94820f4b53cSBarry Smith   relating the cells of the pre-adaptation forest to the post-adaptiation forest.
949bf9b5d84SToby Isaac 
95020f4b53cSBarry Smith   Logically Collective
951bf9b5d84SToby Isaac 
952bf9b5d84SToby Isaac   Input Parameters:
953bf9b5d84SToby Isaac + dm        - the post-adaptation forest
95420f4b53cSBarry Smith - computeSF - default `PETSC_TRUE`
955bf9b5d84SToby Isaac 
956bf9b5d84SToby Isaac   Level: advanced
957bf9b5d84SToby Isaac 
95820f4b53cSBarry Smith   Note:
95920f4b53cSBarry Smith   After `DMSetUp()` is called, the transfer `PetscSF`s can be accessed with `DMForestGetAdaptivitySF()`.
96020f4b53cSBarry Smith 
96120f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetComputeAdaptivitySF()`, `DMForestGetAdaptivitySF()`
962bf9b5d84SToby Isaac @*/
963d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
964d71ae5a4SJacob Faibussowitsch {
965bf9b5d84SToby Isaac   DM_Forest *forest;
966bf9b5d84SToby Isaac 
967bf9b5d84SToby Isaac   PetscFunctionBegin;
968bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
96928b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot compute adaptivity PetscSFs after setup is called");
970bf9b5d84SToby Isaac   forest                 = (DM_Forest *)dm->data;
971bf9b5d84SToby Isaac   forest->computeAdaptSF = computeSF;
9723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
973bf9b5d84SToby Isaac }
974bf9b5d84SToby Isaac 
975d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time)
976d71ae5a4SJacob Faibussowitsch {
97780b27e07SToby Isaac   DM_Forest *forest;
97880b27e07SToby Isaac 
97980b27e07SToby Isaac   PetscFunctionBegin;
98080b27e07SToby Isaac   PetscValidHeaderSpecific(dmIn, DM_CLASSID, 1);
98180b27e07SToby Isaac   PetscValidHeaderSpecific(vecIn, VEC_CLASSID, 2);
98280b27e07SToby Isaac   PetscValidHeaderSpecific(dmOut, DM_CLASSID, 3);
98380b27e07SToby Isaac   PetscValidHeaderSpecific(vecOut, VEC_CLASSID, 4);
98480b27e07SToby Isaac   forest = (DM_Forest *)dmIn->data;
98528b400f6SJacob Faibussowitsch   PetscCheck(forest->transfervec, PetscObjectComm((PetscObject)dmIn), PETSC_ERR_SUP, "DMForestTransferVec() not implemented");
9869566063dSJacob Faibussowitsch   PetscCall((forest->transfervec)(dmIn, vecIn, dmOut, vecOut, useBCs, time));
9873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
98880b27e07SToby Isaac }
98980b27e07SToby Isaac 
990d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestTransferVecFromBase(DM dm, Vec vecIn, Vec vecOut)
991d71ae5a4SJacob Faibussowitsch {
992ac34a06fSStefano Zampini   DM_Forest *forest;
993ac34a06fSStefano Zampini 
994ac34a06fSStefano Zampini   PetscFunctionBegin;
995ac34a06fSStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
996ac34a06fSStefano Zampini   PetscValidHeaderSpecific(vecIn, VEC_CLASSID, 2);
997ac34a06fSStefano Zampini   PetscValidHeaderSpecific(vecOut, VEC_CLASSID, 3);
998ac34a06fSStefano Zampini   forest = (DM_Forest *)dm->data;
99928b400f6SJacob Faibussowitsch   PetscCheck(forest->transfervecfrombase, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMForestTransferVecFromBase() not implemented");
10009566063dSJacob Faibussowitsch   PetscCall((forest->transfervecfrombase)(dm, vecIn, vecOut));
10013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1002ac34a06fSStefano Zampini }
1003ac34a06fSStefano Zampini 
1004bf9b5d84SToby Isaac /*@
100520f4b53cSBarry Smith   DMForestGetComputeAdaptivitySF - Get whether transfer `PetscSF`s should be computed relating the cells of the
100620f4b53cSBarry Smith   pre-adaptation forest to the post-adaptiation forest.  After `DMSetUp()` is called, these transfer PetscSFs can be
100720f4b53cSBarry Smith   accessed with `DMForestGetAdaptivitySF()`.
1008bf9b5d84SToby Isaac 
100920f4b53cSBarry Smith   Not Collective
1010bf9b5d84SToby Isaac 
1011bf9b5d84SToby Isaac   Input Parameter:
1012bf9b5d84SToby Isaac . dm - the post-adaptation forest
1013bf9b5d84SToby Isaac 
1014bf9b5d84SToby Isaac   Output Parameter:
101520f4b53cSBarry Smith . computeSF - default `PETSC_TRUE`
1016bf9b5d84SToby Isaac 
1017bf9b5d84SToby Isaac   Level: advanced
1018bf9b5d84SToby Isaac 
101920f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetComputeAdaptivitySF()`, `DMForestGetAdaptivitySF()`
1020bf9b5d84SToby Isaac @*/
1021d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1022d71ae5a4SJacob Faibussowitsch {
1023bf9b5d84SToby Isaac   DM_Forest *forest;
1024bf9b5d84SToby Isaac 
1025bf9b5d84SToby Isaac   PetscFunctionBegin;
1026bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1027bf9b5d84SToby Isaac   forest     = (DM_Forest *)dm->data;
1028bf9b5d84SToby Isaac   *computeSF = forest->computeAdaptSF;
10293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1030bf9b5d84SToby Isaac }
1031bf9b5d84SToby Isaac 
1032bf9b5d84SToby Isaac /*@
103320f4b53cSBarry Smith   DMForestGetAdaptivitySF - Get `PetscSF`s that relate the pre-adaptation forest to the post-adaptation forest.
1034bf9b5d84SToby Isaac   Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be
1035bf9b5d84SToby Isaac   some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa.  Therefore there are two
103620f4b53cSBarry Smith   `PetscSF`s: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates
1037bf9b5d84SToby Isaac   pre-adaptation fine cells to post-adaptation coarse cells.
1038bf9b5d84SToby Isaac 
103920f4b53cSBarry Smith   Not Collective
1040bf9b5d84SToby Isaac 
1041bf9b5d84SToby Isaac   Input Parameter:
104220f4b53cSBarry Smith . dm - the post-adaptation forest
1043bf9b5d84SToby Isaac 
104420f4b53cSBarry Smith   Output Parameters:
104520f4b53cSBarry Smith + preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
104620f4b53cSBarry Smith - coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-
1047bf9b5d84SToby Isaac 
1048bf9b5d84SToby Isaac   Level: advanced
1049bf9b5d84SToby Isaac 
105020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetComputeAdaptivitySF()`, `DMForestSetComputeAdaptivitySF()`
1051bf9b5d84SToby Isaac @*/
1052d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1053d71ae5a4SJacob Faibussowitsch {
1054bf9b5d84SToby Isaac   DM_Forest *forest;
1055bf9b5d84SToby Isaac 
1056bf9b5d84SToby Isaac   PetscFunctionBegin;
1057bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10589566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
1059bf9b5d84SToby Isaac   forest = (DM_Forest *)dm->data;
1060f885a11aSToby Isaac   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1061f885a11aSToby Isaac   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
10623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1063bf9b5d84SToby Isaac }
1064bf9b5d84SToby Isaac 
10659be51f97SToby Isaac /*@
10669be51f97SToby Isaac   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
106720f4b53cSBarry Smith   indicate that the diameter of neighboring cells should differ by at most a factor of 2.  Subtypes of `DMFOREST` may
10689be51f97SToby Isaac   only support one particular choice of grading factor.
10699be51f97SToby Isaac 
107020f4b53cSBarry Smith   Logically Collective
10719be51f97SToby Isaac 
10729be51f97SToby Isaac   Input Parameters:
10739be51f97SToby Isaac + dm    - the forest
10749be51f97SToby Isaac - grade - the grading factor
10759be51f97SToby Isaac 
10769be51f97SToby Isaac   Level: advanced
10779be51f97SToby Isaac 
107820f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetGradeFactor()`
10799be51f97SToby Isaac @*/
1080d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1081d71ae5a4SJacob Faibussowitsch {
1082c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1083c7eeac06SToby Isaac 
1084c7eeac06SToby Isaac   PetscFunctionBegin;
1085c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
108628b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the grade factor after setup");
1087c7eeac06SToby Isaac   forest->gradeFactor = grade;
10883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1089c7eeac06SToby Isaac }
1090c7eeac06SToby Isaac 
10919be51f97SToby Isaac /*@
10929be51f97SToby Isaac   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
109320f4b53cSBarry Smith   neighboring cells should differ by at most a factor of 2.  Subtypes of `DMFOREST` may only support one particular
10949be51f97SToby Isaac   choice of grading factor.
10959be51f97SToby Isaac 
109620f4b53cSBarry Smith   Not Collective
10979be51f97SToby Isaac 
10989be51f97SToby Isaac   Input Parameter:
10999be51f97SToby Isaac . dm - the forest
11009be51f97SToby Isaac 
11019be51f97SToby Isaac   Output Parameter:
11029be51f97SToby Isaac . grade - the grading factor
11039be51f97SToby Isaac 
11049be51f97SToby Isaac   Level: advanced
11059be51f97SToby Isaac 
110620f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetGradeFactor()`
11079be51f97SToby Isaac @*/
1108d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1109d71ae5a4SJacob Faibussowitsch {
1110c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1111c7eeac06SToby Isaac 
1112c7eeac06SToby Isaac   PetscFunctionBegin;
1113c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1114c7eeac06SToby Isaac   PetscValidIntPointer(grade, 2);
1115c7eeac06SToby Isaac   *grade = forest->gradeFactor;
11163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1117c7eeac06SToby Isaac }
1118c7eeac06SToby Isaac 
11199be51f97SToby Isaac /*@
11209be51f97SToby Isaac   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
112120f4b53cSBarry Smith   the cell weight (see `DMForestSetCellWeights()`) when calculating partitions.  The final weight of a cell will be
11229be51f97SToby Isaac   (cellWeight) * (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on
11239be51f97SToby Isaac   its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
11249be51f97SToby Isaac   computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11259be51f97SToby Isaac 
112620f4b53cSBarry Smith   Logically Collective
11279be51f97SToby Isaac 
11289be51f97SToby Isaac   Input Parameters:
11299be51f97SToby Isaac + dm            - the forest
1130*60225df5SJacob Faibussowitsch - weightsFactor - default 1.
11319be51f97SToby Isaac 
11329be51f97SToby Isaac   Level: advanced
11339be51f97SToby Isaac 
113420f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetCellWeightFactor()`, `DMForestSetCellWeights()`
11359be51f97SToby Isaac @*/
1136d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
1137d71ae5a4SJacob Faibussowitsch {
1138c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1139c7eeac06SToby Isaac 
1140c7eeac06SToby Isaac   PetscFunctionBegin;
1141c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
114228b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the weights factor after setup");
1143c7eeac06SToby Isaac   forest->weightsFactor = weightsFactor;
11443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1145c7eeac06SToby Isaac }
1146c7eeac06SToby Isaac 
11479be51f97SToby Isaac /*@
11489be51f97SToby Isaac   DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
114920f4b53cSBarry Smith   `DMForestSetCellWeights()`) when calculating partitions.  The final weight of a cell will be (cellWeight) *
11509be51f97SToby Isaac   (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on its level; a
11519be51f97SToby Isaac   factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation
11529be51f97SToby Isaac   associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11539be51f97SToby Isaac 
115420f4b53cSBarry Smith   Not Collective
11559be51f97SToby Isaac 
11569be51f97SToby Isaac   Input Parameter:
11579be51f97SToby Isaac . dm - the forest
11589be51f97SToby Isaac 
11599be51f97SToby Isaac   Output Parameter:
1160*60225df5SJacob Faibussowitsch . weightsFactor - default 1.
11619be51f97SToby Isaac 
11629be51f97SToby Isaac   Level: advanced
11639be51f97SToby Isaac 
116420f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetCellWeightFactor()`, `DMForestSetCellWeights()`
11659be51f97SToby Isaac @*/
1166d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1167d71ae5a4SJacob Faibussowitsch {
1168c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1169c7eeac06SToby Isaac 
1170c7eeac06SToby Isaac   PetscFunctionBegin;
1171c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1172c7eeac06SToby Isaac   PetscValidRealPointer(weightsFactor, 2);
1173c7eeac06SToby Isaac   *weightsFactor = forest->weightsFactor;
11743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1175c7eeac06SToby Isaac }
1176c7eeac06SToby Isaac 
11779be51f97SToby Isaac /*@
11789be51f97SToby Isaac   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
11799be51f97SToby Isaac 
118020f4b53cSBarry Smith   Not Collective
11819be51f97SToby Isaac 
11829be51f97SToby Isaac   Input Parameter:
11839be51f97SToby Isaac . dm - the forest
11849be51f97SToby Isaac 
11859be51f97SToby Isaac   Output Parameters:
11869be51f97SToby Isaac + cStart - the first cell on this process
11879be51f97SToby Isaac - cEnd   - one after the final cell on this process
11889be51f97SToby Isaac 
11891a244344SSatish Balay   Level: intermediate
11909be51f97SToby Isaac 
119120f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetCellSF()`
11929be51f97SToby Isaac @*/
1193d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1194d71ae5a4SJacob Faibussowitsch {
1195c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1196c7eeac06SToby Isaac 
1197c7eeac06SToby Isaac   PetscFunctionBegin;
1198c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1199c7eeac06SToby Isaac   PetscValidIntPointer(cStart, 2);
1200064a246eSJacob Faibussowitsch   PetscValidIntPointer(cEnd, 3);
120148a46eb9SPierre Jolivet   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) PetscCall(forest->createcellchart(dm, &forest->cStart, &forest->cEnd));
1202c7eeac06SToby Isaac   *cStart = forest->cStart;
1203c7eeac06SToby Isaac   *cEnd   = forest->cEnd;
12043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1205c7eeac06SToby Isaac }
1206c7eeac06SToby Isaac 
12079be51f97SToby Isaac /*@
120820f4b53cSBarry Smith   DMForestGetCellSF - After the setup phase, get the `PetscSF` for overlapping cells between processes
12099be51f97SToby Isaac 
121020f4b53cSBarry Smith   Not Collective
12119be51f97SToby Isaac 
12129be51f97SToby Isaac   Input Parameter:
12139be51f97SToby Isaac . dm - the forest
12149be51f97SToby Isaac 
12159be51f97SToby Isaac   Output Parameter:
121620f4b53cSBarry Smith . cellSF - the `PetscSF`
12179be51f97SToby Isaac 
12181a244344SSatish Balay   Level: intermediate
12199be51f97SToby Isaac 
122020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetCellChart()`
12219be51f97SToby Isaac @*/
1222d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1223d71ae5a4SJacob Faibussowitsch {
1224c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1225c7eeac06SToby Isaac 
1226c7eeac06SToby Isaac   PetscFunctionBegin;
1227c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1228c7eeac06SToby Isaac   PetscValidPointer(cellSF, 2);
122948a46eb9SPierre Jolivet   if ((!forest->cellSF) && forest->createcellsf) PetscCall(forest->createcellsf(dm, &forest->cellSF));
1230c7eeac06SToby Isaac   *cellSF = forest->cellSF;
12313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1232c7eeac06SToby Isaac }
1233c7eeac06SToby Isaac 
12349be51f97SToby Isaac /*@C
12359be51f97SToby Isaac   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
123620f4b53cSBarry Smith   `DMForestGetAdaptivityForest()`) that holds the adaptation flags (refinement, coarsening, or some combination).  The
123720f4b53cSBarry Smith   interpretation of the label values is up to the subtype of `DMFOREST`, but `DM_ADAPT_DETERMINE`, `DM_ADAPT_KEEP`,
123820f4b53cSBarry Smith   `DM_ADAPT_REFINE`, and `DM_ADAPT_COARSEN` have been reserved as choices that should be accepted by all subtypes.
12399be51f97SToby Isaac 
124020f4b53cSBarry Smith   Logically Collective
12419be51f97SToby Isaac 
12429be51f97SToby Isaac   Input Parameters:
1243*60225df5SJacob Faibussowitsch + dm         - the forest
1244*60225df5SJacob Faibussowitsch - adaptLabel - the label in the pre-adaptation forest
12459be51f97SToby Isaac 
12469be51f97SToby Isaac   Level: intermediate
12479be51f97SToby Isaac 
124820f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdaptivityLabel()`
12499be51f97SToby Isaac @*/
1250d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel)
1251d71ae5a4SJacob Faibussowitsch {
1252c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1253c7eeac06SToby Isaac 
1254c7eeac06SToby Isaac   PetscFunctionBegin;
1255c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12561fd50544SStefano Zampini   if (adaptLabel) PetscValidHeaderSpecific(adaptLabel, DMLABEL_CLASSID, 2);
12579566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)adaptLabel));
12589566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&forest->adaptLabel));
12591fd50544SStefano Zampini   forest->adaptLabel = adaptLabel;
12603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1261c7eeac06SToby Isaac }
1262c7eeac06SToby Isaac 
12639be51f97SToby Isaac /*@C
126420f4b53cSBarry Smith   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see `DMForestGetAdaptivityForest()`) that
12659be51f97SToby Isaac   holds the adaptation flags (refinement, coarsening, or some combination).  The interpretation of the label values is
126620f4b53cSBarry Smith   up to the subtype of `DMFOREST`, but `DM_ADAPT_DETERMINE`, `DM_ADAPT_KEEP`, `DM_ADAPT_REFINE`, and `DM_ADAPT_COARSEN` have
1267cd3c525cSToby Isaac   been reserved as choices that should be accepted by all subtypes.
12689be51f97SToby Isaac 
126920f4b53cSBarry Smith   Not Collective
12709be51f97SToby Isaac 
12719be51f97SToby Isaac   Input Parameter:
12729be51f97SToby Isaac . dm - the forest
12739be51f97SToby Isaac 
12749be51f97SToby Isaac   Output Parameter:
12759be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest
12769be51f97SToby Isaac 
12779be51f97SToby Isaac   Level: intermediate
12789be51f97SToby Isaac 
127920f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdaptivityLabel()`
12809be51f97SToby Isaac @*/
1281d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel)
1282d71ae5a4SJacob Faibussowitsch {
1283c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1284c7eeac06SToby Isaac 
1285c7eeac06SToby Isaac   PetscFunctionBegin;
1286c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1287ba936b91SToby Isaac   *adaptLabel = forest->adaptLabel;
12883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1289c7eeac06SToby Isaac }
1290c7eeac06SToby Isaac 
12919be51f97SToby Isaac /*@
129220f4b53cSBarry Smith   DMForestSetCellWeights - Set the weights assigned to each of the cells (see `DMForestGetCellChart()`) of the current
129320f4b53cSBarry Smith   process: weights are used to determine parallel partitioning.
12949be51f97SToby Isaac 
129520f4b53cSBarry Smith   Logically Collective
12969be51f97SToby Isaac 
12979be51f97SToby Isaac   Input Parameters:
12989be51f97SToby Isaac + dm       - the forest
129920f4b53cSBarry Smith . weights  - the array of weights (see `DMForestSetWeightCapacity()`) for all cells, or `NULL` to indicate each cell has weight 1.
13009be51f97SToby Isaac - copyMode - how weights should reference weights
13019be51f97SToby Isaac 
13029be51f97SToby Isaac   Level: advanced
13039be51f97SToby Isaac 
130420f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetCellWeights()`, `DMForestSetWeightCapacity()`
13059be51f97SToby Isaac @*/
1306d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1307d71ae5a4SJacob Faibussowitsch {
1308c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1309c7eeac06SToby Isaac   PetscInt   cStart, cEnd;
1310c7eeac06SToby Isaac 
1311c7eeac06SToby Isaac   PetscFunctionBegin;
1312c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
13139566063dSJacob Faibussowitsch   PetscCall(DMForestGetCellChart(dm, &cStart, &cEnd));
131463a3b9bcSJacob Faibussowitsch   PetscCheck(cEnd >= cStart, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "cell chart [%" PetscInt_FMT ",%" PetscInt_FMT ") is not valid", cStart, cEnd);
1315c7eeac06SToby Isaac   if (copyMode == PETSC_COPY_VALUES) {
131648a46eb9SPierre Jolivet     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) PetscCall(PetscMalloc1(cEnd - cStart, &forest->cellWeights));
13179566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(forest->cellWeights, weights, cEnd - cStart));
1318c7eeac06SToby Isaac     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
13193ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1320c7eeac06SToby Isaac   }
132148a46eb9SPierre Jolivet   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) PetscCall(PetscFree(forest->cellWeights));
1322c7eeac06SToby Isaac   forest->cellWeights         = weights;
1323c7eeac06SToby Isaac   forest->cellWeightsCopyMode = copyMode;
13243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1325c7eeac06SToby Isaac }
1326c7eeac06SToby Isaac 
13279be51f97SToby Isaac /*@
132820f4b53cSBarry Smith   DMForestGetCellWeights - Get the weights assigned to each of the cells (see `DMForestGetCellChart()`) of the current
132920f4b53cSBarry Smith   process: weights are used to determine parallel partitioning.
13309be51f97SToby Isaac 
133120f4b53cSBarry Smith   Not Collective
13329be51f97SToby Isaac 
13339be51f97SToby Isaac   Input Parameter:
13349be51f97SToby Isaac . dm - the forest
13359be51f97SToby Isaac 
13369be51f97SToby Isaac   Output Parameter:
133720f4b53cSBarry Smith . weights - the array of weights for all cells, or `NULL` to indicate each cell has weight 1.
13389be51f97SToby Isaac 
13399be51f97SToby Isaac   Level: advanced
13409be51f97SToby Isaac 
134120f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetCellWeights()`, `DMForestSetWeightCapacity()`
13429be51f97SToby Isaac @*/
1343d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1344d71ae5a4SJacob Faibussowitsch {
1345c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1346c7eeac06SToby Isaac 
1347c7eeac06SToby Isaac   PetscFunctionBegin;
1348c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1349c7eeac06SToby Isaac   PetscValidPointer(weights, 2);
1350c7eeac06SToby Isaac   *weights = forest->cellWeights;
13513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1352c7eeac06SToby Isaac }
1353c7eeac06SToby Isaac 
13549be51f97SToby Isaac /*@
13559be51f97SToby Isaac   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
135620f4b53cSBarry Smith   a pre-adaptation forest (see `DMForestGetAdaptivityForest()`).  After partitioning, the ratio of the weight of each
13579be51f97SToby Isaac   process's cells to the process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that
13589be51f97SToby Isaac   the current process should not have any cells after repartitioning.
13599be51f97SToby Isaac 
136020f4b53cSBarry Smith   Logically Collective
13619be51f97SToby Isaac 
1362*60225df5SJacob Faibussowitsch   Input Parameters:
13639be51f97SToby Isaac + dm       - the forest
13649be51f97SToby Isaac - capacity - this process's capacity
13659be51f97SToby Isaac 
13669be51f97SToby Isaac   Level: advanced
13679be51f97SToby Isaac 
136820f4b53cSBarry Smith .seealso `DM`, `DMFOREST`, `DMForestGetWeightCapacity()`, `DMForestSetCellWeights()`, `DMForestSetCellWeightFactor()`
13699be51f97SToby Isaac @*/
1370d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1371d71ae5a4SJacob Faibussowitsch {
1372c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1373c7eeac06SToby Isaac 
1374c7eeac06SToby Isaac   PetscFunctionBegin;
1375c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
137628b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the weight capacity after setup");
137763a3b9bcSJacob Faibussowitsch   PetscCheck(capacity >= 0., PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot have negative weight capacity; %g", (double)capacity);
1378c7eeac06SToby Isaac   forest->weightCapacity = capacity;
13793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1380c7eeac06SToby Isaac }
1381c7eeac06SToby Isaac 
13829be51f97SToby Isaac /*@
13839be51f97SToby Isaac   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
138420f4b53cSBarry Smith   `DMForestGetAdaptivityForest()`).  After partitioning, the ratio of the weight of each process's cells to the
13859be51f97SToby Isaac   process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that the current process
13869be51f97SToby Isaac   should not have any cells after repartitioning.
13879be51f97SToby Isaac 
138820f4b53cSBarry Smith   Not Collective
13899be51f97SToby Isaac 
1390*60225df5SJacob Faibussowitsch   Input Parameter:
13919be51f97SToby Isaac . dm - the forest
13929be51f97SToby Isaac 
1393*60225df5SJacob Faibussowitsch   Output Parameter:
13949be51f97SToby Isaac . capacity - this process's capacity
13959be51f97SToby Isaac 
13969be51f97SToby Isaac   Level: advanced
13979be51f97SToby Isaac 
139820f4b53cSBarry Smith .seealso `DM`, `DMFOREST`, `DMForestSetWeightCapacity()`, `DMForestSetCellWeights()`, `DMForestSetCellWeightFactor()`
13999be51f97SToby Isaac @*/
1400d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1401d71ae5a4SJacob Faibussowitsch {
1402c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1403c7eeac06SToby Isaac 
1404c7eeac06SToby Isaac   PetscFunctionBegin;
1405c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1406c7eeac06SToby Isaac   PetscValidRealPointer(capacity, 2);
1407c7eeac06SToby Isaac   *capacity = forest->weightCapacity;
14083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1409c7eeac06SToby Isaac }
1410c7eeac06SToby Isaac 
1411d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(DM dm, PetscOptionItems *PetscOptionsObject)
1412d71ae5a4SJacob Faibussowitsch {
141356ba9f64SToby Isaac   PetscBool                  flg, flg1, flg2, flg3, flg4;
1414dd8e54a2SToby Isaac   DMForestTopology           oldTopo;
1415c7eeac06SToby Isaac   char                       stringBuffer[256];
1416dd8e54a2SToby Isaac   PetscViewer                viewer;
1417dd8e54a2SToby Isaac   PetscViewerFormat          format;
141856ba9f64SToby Isaac   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1419c7eeac06SToby Isaac   PetscReal                  weightsFactor;
1420c7eeac06SToby Isaac   DMForestAdaptivityStrategy adaptStrategy;
1421db4d5e8cSToby Isaac 
1422db4d5e8cSToby Isaac   PetscFunctionBegin;
14239566063dSJacob Faibussowitsch   PetscCall(DMForestGetTopology(dm, &oldTopo));
1424d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "DMForest Options");
14259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-dm_forest_topology", "the topology of the forest's base mesh", "DMForestSetTopology", oldTopo, stringBuffer, sizeof(stringBuffer), &flg1));
14269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsViewer("-dm_forest_base_dm", "load the base DM from a viewer specification", "DMForestSetBaseDM", &viewer, &format, &flg2));
14279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsViewer("-dm_forest_coarse_forest", "load the coarse forest from a viewer specification", "DMForestSetCoarseForest", &viewer, &format, &flg3));
14289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsViewer("-dm_forest_fine_forest", "load the fine forest from a viewer specification", "DMForestSetFineForest", &viewer, &format, &flg4));
14291dca8a05SBarry 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}");
143056ba9f64SToby Isaac   if (flg1) {
14319566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm, (DMForestTopology)stringBuffer));
14329566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm, NULL));
14339566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, NULL));
143456ba9f64SToby Isaac   }
143556ba9f64SToby Isaac   if (flg2) {
1436dd8e54a2SToby Isaac     DM base;
1437dd8e54a2SToby Isaac 
14389566063dSJacob Faibussowitsch     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &base));
14399566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
14409566063dSJacob Faibussowitsch     PetscCall(DMLoad(base, viewer));
14419566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
14429566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm, base));
14439566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&base));
14449566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm, NULL));
14459566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, NULL));
1446dd8e54a2SToby Isaac   }
144756ba9f64SToby Isaac   if (flg3) {
1448dd8e54a2SToby Isaac     DM coarse;
1449dd8e54a2SToby Isaac 
14509566063dSJacob Faibussowitsch     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &coarse));
14519566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
14529566063dSJacob Faibussowitsch     PetscCall(DMLoad(coarse, viewer));
14539566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
14549566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, coarse));
14559566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&coarse));
14569566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm, NULL));
14579566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm, NULL));
1458dd8e54a2SToby Isaac   }
145956ba9f64SToby Isaac   if (flg4) {
1460dd8e54a2SToby Isaac     DM fine;
1461dd8e54a2SToby Isaac 
14629566063dSJacob Faibussowitsch     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &fine));
14639566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
14649566063dSJacob Faibussowitsch     PetscCall(DMLoad(fine, viewer));
14659566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
14669566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, fine));
14679566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&fine));
14689566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm, NULL));
14699566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm, NULL));
1470dd8e54a2SToby Isaac   }
14719566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdjacencyDimension(dm, &adjDim));
14729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_dimension", "set the dimension of points that define adjacency in the forest", "DMForestSetAdjacencyDimension", adjDim, &adjDim, &flg, 0));
1473dd8e54a2SToby Isaac   if (flg) {
14749566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdjacencyDimension(dm, adjDim));
1475f885a11aSToby Isaac   } else {
14769566063dSJacob Faibussowitsch     PetscCall(DMForestGetAdjacencyCodimension(dm, &adjCodim));
14779566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_codimension", "set the codimension of points that define adjacency in the forest", "DMForestSetAdjacencyCodimension", adjCodim, &adjCodim, &flg, 1));
14781baa6e33SBarry Smith     if (flg) PetscCall(DMForestSetAdjacencyCodimension(dm, adjCodim));
1479dd8e54a2SToby Isaac   }
14809566063dSJacob Faibussowitsch   PetscCall(DMForestGetPartitionOverlap(dm, &overlap));
14819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_partition_overlap", "set the degree of partition overlap", "DMForestSetPartitionOverlap", overlap, &overlap, &flg, 0));
14821baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetPartitionOverlap(dm, overlap));
1483a6121fbdSMatthew G. Knepley #if 0
14849566063dSJacob 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));
1485a6121fbdSMatthew G. Knepley   if (flg) {
14869566063dSJacob Faibussowitsch     PetscCall(DMForestSetMinimumRefinement(dm,minRefinement));
14879566063dSJacob Faibussowitsch     PetscCall(DMForestSetInitialRefinement(dm,minRefinement));
1488a6121fbdSMatthew G. Knepley   }
14899566063dSJacob 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));
1490a6121fbdSMatthew G. Knepley   if (flg) {
14919566063dSJacob Faibussowitsch     PetscCall(DMForestSetMinimumRefinement(dm,0));
14929566063dSJacob Faibussowitsch     PetscCall(DMForestSetInitialRefinement(dm,initRefinement));
1493a6121fbdSMatthew G. Knepley   }
1494a6121fbdSMatthew G. Knepley #endif
14959566063dSJacob Faibussowitsch   PetscCall(DMForestGetMinimumRefinement(dm, &minRefinement));
14969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_minimum_refinement", "set the minimum level of refinement in the forest", "DMForestSetMinimumRefinement", minRefinement, &minRefinement, &flg, 0));
14971baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetMinimumRefinement(dm, minRefinement));
14989566063dSJacob Faibussowitsch   PetscCall(DMForestGetInitialRefinement(dm, &initRefinement));
14999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_initial_refinement", "set the initial level of refinement in the forest", "DMForestSetInitialRefinement", initRefinement, &initRefinement, &flg, 0));
15001baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetInitialRefinement(dm, initRefinement));
15019566063dSJacob Faibussowitsch   PetscCall(DMForestGetMaximumRefinement(dm, &maxRefinement));
15029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_maximum_refinement", "set the maximum level of refinement in the forest", "DMForestSetMaximumRefinement", maxRefinement, &maxRefinement, &flg, 0));
15031baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetMaximumRefinement(dm, maxRefinement));
15049566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivityStrategy(dm, &adaptStrategy));
15059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-dm_forest_adaptivity_strategy", "the forest's adaptivity-flag resolution strategy", "DMForestSetAdaptivityStrategy", adaptStrategy, stringBuffer, sizeof(stringBuffer), &flg));
15061baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetAdaptivityStrategy(dm, (DMForestAdaptivityStrategy)stringBuffer));
15079566063dSJacob Faibussowitsch   PetscCall(DMForestGetGradeFactor(dm, &grade));
15089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_grade_factor", "grade factor between neighboring cells", "DMForestSetGradeFactor", grade, &grade, &flg, 0));
15091baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetGradeFactor(dm, grade));
15109566063dSJacob Faibussowitsch   PetscCall(DMForestGetCellWeightFactor(dm, &weightsFactor));
15119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_forest_cell_weight_factor", "multiplying weight factor for cell refinement", "DMForestSetCellWeightFactor", weightsFactor, &weightsFactor, &flg));
15121baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetCellWeightFactor(dm, weightsFactor));
1513d0609cedSBarry Smith   PetscOptionsHeadEnd();
15143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1515db4d5e8cSToby Isaac }
1516db4d5e8cSToby Isaac 
1517d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1518d71ae5a4SJacob Faibussowitsch {
1519d8984e3bSMatthew G. Knepley   PetscFunctionBegin;
15209566063dSJacob Faibussowitsch   if (subdm) PetscCall(DMClone(dm, subdm));
15219566063dSJacob Faibussowitsch   PetscCall(DMCreateSectionSubDM(dm, numFields, fields, is, subdm));
15223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1523d8984e3bSMatthew G. Knepley }
1524d8984e3bSMatthew G. Knepley 
1525d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
1526d71ae5a4SJacob Faibussowitsch {
15275421bac9SToby Isaac   DMLabel refine;
15285421bac9SToby Isaac   DM      fineDM;
15295421bac9SToby Isaac 
15305421bac9SToby Isaac   PetscFunctionBegin;
15319566063dSJacob Faibussowitsch   PetscCall(DMGetFineDM(dm, &fineDM));
15325421bac9SToby Isaac   if (fineDM) {
15339566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fineDM));
15345421bac9SToby Isaac     *dmRefined = fineDM;
15353ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
15365421bac9SToby Isaac   }
15379566063dSJacob Faibussowitsch   PetscCall(DMForestTemplate(dm, comm, dmRefined));
15389566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "refine", &refine));
15395421bac9SToby Isaac   if (!refine) {
15409566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "refine", &refine));
15419566063dSJacob Faibussowitsch     PetscCall(DMLabelSetDefaultValue(refine, DM_ADAPT_REFINE));
15421baa6e33SBarry Smith   } else PetscCall(PetscObjectReference((PetscObject)refine));
15439566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityLabel(*dmRefined, refine));
15449566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&refine));
15453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15465421bac9SToby Isaac }
15475421bac9SToby Isaac 
1548d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
1549d71ae5a4SJacob Faibussowitsch {
15505421bac9SToby Isaac   DMLabel coarsen;
15515421bac9SToby Isaac   DM      coarseDM;
15525421bac9SToby Isaac 
15535421bac9SToby Isaac   PetscFunctionBegin;
15544098eed7SToby Isaac   {
15554098eed7SToby Isaac     PetscMPIInt mpiComparison;
15564098eed7SToby Isaac     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);
15574098eed7SToby Isaac 
15589566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_compare(comm, dmcomm, &mpiComparison));
15591dca8a05SBarry Smith     PetscCheck(mpiComparison == MPI_IDENT || mpiComparison == MPI_CONGRUENT, dmcomm, PETSC_ERR_SUP, "No support for different communicators yet");
15604098eed7SToby Isaac   }
15619566063dSJacob Faibussowitsch   PetscCall(DMGetCoarseDM(dm, &coarseDM));
15625421bac9SToby Isaac   if (coarseDM) {
15639566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)coarseDM));
15645421bac9SToby Isaac     *dmCoarsened = coarseDM;
15653ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
15665421bac9SToby Isaac   }
15679566063dSJacob Faibussowitsch   PetscCall(DMForestTemplate(dm, comm, dmCoarsened));
15689566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityPurpose(*dmCoarsened, DM_ADAPT_COARSEN));
15699566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "coarsen", &coarsen));
15705421bac9SToby Isaac   if (!coarsen) {
15719566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "coarsen", &coarsen));
15729566063dSJacob Faibussowitsch     PetscCall(DMLabelSetDefaultValue(coarsen, DM_ADAPT_COARSEN));
15731baa6e33SBarry Smith   } else PetscCall(PetscObjectReference((PetscObject)coarsen));
15749566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityLabel(*dmCoarsened, coarsen));
15759566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&coarsen));
15763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15775421bac9SToby Isaac }
15785421bac9SToby Isaac 
1579d71ae5a4SJacob Faibussowitsch PetscErrorCode DMAdaptLabel_Forest(DM dm, PETSC_UNUSED Vec metric, DMLabel label, PETSC_UNUSED DMLabel rgLabel, DM *adaptedDM)
1580d71ae5a4SJacob Faibussowitsch {
158109350103SToby Isaac   PetscBool success;
158209350103SToby Isaac 
158309350103SToby Isaac   PetscFunctionBegin;
15849566063dSJacob Faibussowitsch   PetscCall(DMForestTemplate(dm, PetscObjectComm((PetscObject)dm), adaptedDM));
15859566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityLabel(*adaptedDM, label));
15869566063dSJacob Faibussowitsch   PetscCall(DMSetUp(*adaptedDM));
15879566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivitySuccess(*adaptedDM, &success));
158809350103SToby Isaac   if (!success) {
15899566063dSJacob Faibussowitsch     PetscCall(DMDestroy(adaptedDM));
159009350103SToby Isaac     *adaptedDM = NULL;
159109350103SToby Isaac   }
15923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
159309350103SToby Isaac }
159409350103SToby Isaac 
1595d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMInitialize_Forest(DM dm)
1596d71ae5a4SJacob Faibussowitsch {
1597d222f98bSToby Isaac   PetscFunctionBegin;
15989566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(dm->ops, sizeof(*(dm->ops))));
1599d222f98bSToby Isaac 
1600d222f98bSToby Isaac   dm->ops->clone          = DMClone_Forest;
1601d222f98bSToby Isaac   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1602d222f98bSToby Isaac   dm->ops->destroy        = DMDestroy_Forest;
1603d8984e3bSMatthew G. Knepley   dm->ops->createsubdm    = DMCreateSubDM_Forest;
16045421bac9SToby Isaac   dm->ops->refine         = DMRefine_Forest;
16055421bac9SToby Isaac   dm->ops->coarsen        = DMCoarsen_Forest;
16063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1607d222f98bSToby Isaac }
1608d222f98bSToby Isaac 
16099be51f97SToby Isaac /*MC
16109be51f97SToby Isaac 
161120f4b53cSBarry Smith      DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh.  Forests usually have a base `DM`
161220f4b53cSBarry Smith   (see `DMForestGetBaseDM()`), from which it is refined.  The refinement and partitioning of forests is considered
161320f4b53cSBarry Smith   immutable after `DMSetUp()` is called.  To adapt a mesh, one should call `DMForestTemplate()` to create a new mesh that
161420f4b53cSBarry Smith   will default to being identical to it, specify how that mesh should differ, and then calling `DMSetUp()` on the new
1615bae1f979SBarry Smith   mesh.
1616bae1f979SBarry Smith 
1617bae1f979SBarry Smith   To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the
161820f4b53cSBarry Smith   previous mesh whose values indicate which cells should be refined (`DM_ADAPT_REFINE`) or coarsened (`DM_ADAPT_COARSEN`)
1619bae1f979SBarry Smith   and how (subtypes are free to allow additional values for things like anisotropic refinement).  The label should be
162020f4b53cSBarry Smith   given to the *new* mesh with `DMForestSetAdaptivityLabel()`.
16219be51f97SToby Isaac 
16229be51f97SToby Isaac   Level: advanced
16239be51f97SToby Isaac 
162420f4b53cSBarry Smith .seealso: `DMType`, `DM`, `DMCreate()`, `DMSetType()`, `DMForestGetBaseDM()`, `DMForestSetBaseDM()`, `DMForestTemplate()`, `DMForestSetAdaptivityLabel()`
16259be51f97SToby Isaac M*/
16269be51f97SToby Isaac 
1627d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1628d71ae5a4SJacob Faibussowitsch {
1629db4d5e8cSToby Isaac   DM_Forest *forest;
1630db4d5e8cSToby Isaac 
1631db4d5e8cSToby Isaac   PetscFunctionBegin;
1632db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
16334dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&forest));
1634db4d5e8cSToby Isaac   dm->dim                     = 0;
1635db4d5e8cSToby Isaac   dm->data                    = forest;
1636db4d5e8cSToby Isaac   forest->refct               = 1;
1637db4d5e8cSToby Isaac   forest->data                = NULL;
1638db4d5e8cSToby Isaac   forest->topology            = NULL;
1639cd3c525cSToby Isaac   forest->adapt               = NULL;
1640db4d5e8cSToby Isaac   forest->base                = NULL;
16416a87ffbfSToby Isaac   forest->adaptPurpose        = DM_ADAPT_DETERMINE;
1642db4d5e8cSToby Isaac   forest->adjDim              = PETSC_DEFAULT;
1643db4d5e8cSToby Isaac   forest->overlap             = PETSC_DEFAULT;
1644db4d5e8cSToby Isaac   forest->minRefinement       = PETSC_DEFAULT;
1645db4d5e8cSToby Isaac   forest->maxRefinement       = PETSC_DEFAULT;
164656ba9f64SToby Isaac   forest->initRefinement      = PETSC_DEFAULT;
1647c7eeac06SToby Isaac   forest->cStart              = PETSC_DETERMINE;
1648c7eeac06SToby Isaac   forest->cEnd                = PETSC_DETERMINE;
1649cd3c525cSToby Isaac   forest->cellSF              = NULL;
1650ebdf65a2SToby Isaac   forest->adaptLabel          = NULL;
1651db4d5e8cSToby Isaac   forest->gradeFactor         = 2;
1652db4d5e8cSToby Isaac   forest->cellWeights         = NULL;
1653db4d5e8cSToby Isaac   forest->cellWeightsCopyMode = PETSC_USE_POINTER;
1654db4d5e8cSToby Isaac   forest->weightsFactor       = 1.;
1655db4d5e8cSToby Isaac   forest->weightCapacity      = 1.;
16569566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityStrategy(dm, DMFORESTADAPTALL));
16579566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Forest(dm));
16583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1659db4d5e8cSToby Isaac }
1660