xref: /petsc/src/dm/impls/forest/forest.c (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
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 
4760225df5SJacob 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 
7260225df5SJacob Faibussowitsch   Input Parameter:
739be51f97SToby Isaac . dm - the DM object
749be51f97SToby Isaac 
7560225df5SJacob 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 /*@
101*a4e35b19SJacob Faibussowitsch   DMForestTemplate - Create a new `DM` that will be adapted from a source `DM`.
1029be51f97SToby Isaac 
10320f4b53cSBarry Smith   Collective
1049be51f97SToby Isaac 
1059be51f97SToby Isaac   Input Parameters:
10620f4b53cSBarry Smith + dm   - the source `DM` object
10720f4b53cSBarry Smith - comm - the communicator for the new `DM` (this communicator is currently ignored, but is present so that `DMForestTemplate()` can be used within `DMCoarsen()`)
1089be51f97SToby Isaac 
1099be51f97SToby Isaac   Output Parameter:
11020f4b53cSBarry Smith . tdm - the new `DM` object
1119be51f97SToby Isaac 
1129be51f97SToby Isaac   Level: intermediate
1139be51f97SToby Isaac 
114*a4e35b19SJacob Faibussowitsch   Notes:
115*a4e35b19SJacob Faibussowitsch   The new `DM` reproduces the configuration of the source, but is not yet setup, so that the
116*a4e35b19SJacob Faibussowitsch   user can then define only the ways that the new `DM` should differ (by, e.g., refinement or
117*a4e35b19SJacob Faibussowitsch   repartitioning).  The source `DM` is also set as the adaptivity source `DM` of the new `DM`
118*a4e35b19SJacob Faibussowitsch   (see `DMForestSetAdaptivityForest()`).
119*a4e35b19SJacob Faibussowitsch 
12020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdaptivityForest()`
1219be51f97SToby Isaac @*/
122d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm)
123d71ae5a4SJacob Faibussowitsch {
124a0452a8eSToby Isaac   DM_Forest                 *forest = (DM_Forest *)dm->data;
12520e8089bSToby Isaac   DMType                     type;
126a0452a8eSToby Isaac   DM                         base;
127a0452a8eSToby Isaac   DMForestTopology           topology;
12805e99e11SStefano Zampini   MatType                    mtype;
129a0452a8eSToby Isaac   PetscInt                   dim, overlap, ref, factor;
130a0452a8eSToby Isaac   DMForestAdaptivityStrategy strat;
131795844e7SToby Isaac   void                      *ctx;
13249fc9a2fSToby Isaac   PetscErrorCode (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *);
1333e58adeeSToby Isaac   void *mapCtx;
134a0452a8eSToby Isaac 
135a0452a8eSToby Isaac   PetscFunctionBegin;
136a0452a8eSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1379566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), tdm));
1389566063dSJacob Faibussowitsch   PetscCall(DMGetType(dm, &type));
1399566063dSJacob Faibussowitsch   PetscCall(DMSetType(*tdm, type));
1409566063dSJacob Faibussowitsch   PetscCall(DMForestGetBaseDM(dm, &base));
1419566063dSJacob Faibussowitsch   PetscCall(DMForestSetBaseDM(*tdm, base));
1429566063dSJacob Faibussowitsch   PetscCall(DMForestGetTopology(dm, &topology));
1439566063dSJacob Faibussowitsch   PetscCall(DMForestSetTopology(*tdm, topology));
1449566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdjacencyDimension(dm, &dim));
1459566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdjacencyDimension(*tdm, dim));
1469566063dSJacob Faibussowitsch   PetscCall(DMForestGetPartitionOverlap(dm, &overlap));
1479566063dSJacob Faibussowitsch   PetscCall(DMForestSetPartitionOverlap(*tdm, overlap));
1489566063dSJacob Faibussowitsch   PetscCall(DMForestGetMinimumRefinement(dm, &ref));
1499566063dSJacob Faibussowitsch   PetscCall(DMForestSetMinimumRefinement(*tdm, ref));
1509566063dSJacob Faibussowitsch   PetscCall(DMForestGetMaximumRefinement(dm, &ref));
1519566063dSJacob Faibussowitsch   PetscCall(DMForestSetMaximumRefinement(*tdm, ref));
1529566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivityStrategy(dm, &strat));
1539566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityStrategy(*tdm, strat));
1549566063dSJacob Faibussowitsch   PetscCall(DMForestGetGradeFactor(dm, &factor));
1559566063dSJacob Faibussowitsch   PetscCall(DMForestSetGradeFactor(*tdm, factor));
1569566063dSJacob Faibussowitsch   PetscCall(DMForestGetBaseCoordinateMapping(dm, &map, &mapCtx));
1579566063dSJacob Faibussowitsch   PetscCall(DMForestSetBaseCoordinateMapping(*tdm, map, mapCtx));
1581baa6e33SBarry Smith   if (forest->ftemplate) PetscCall((*forest->ftemplate)(dm, *tdm));
1599566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityForest(*tdm, dm));
1609566063dSJacob Faibussowitsch   PetscCall(DMCopyDisc(dm, *tdm));
1619566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dm, &ctx));
1629566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*tdm, &ctx));
16390b157c4SStefano Zampini   {
1644fb89dddSMatthew G. Knepley     const PetscReal *maxCell, *L, *Lstart;
165795844e7SToby Isaac 
1664fb89dddSMatthew G. Knepley     PetscCall(DMGetPeriodicity(dm, &maxCell, &Lstart, &L));
1674fb89dddSMatthew G. Knepley     PetscCall(DMSetPeriodicity(*tdm, maxCell, Lstart, L));
168795844e7SToby Isaac   }
1699566063dSJacob Faibussowitsch   PetscCall(DMGetMatType(dm, &mtype));
1709566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(*tdm, mtype));
1713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
172a0452a8eSToby Isaac }
173a0452a8eSToby Isaac 
17401d9d024SToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm);
17501d9d024SToby Isaac 
176d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm)
177d71ae5a4SJacob Faibussowitsch {
178db4d5e8cSToby Isaac   DM_Forest  *forest = (DM_Forest *)dm->data;
179db4d5e8cSToby Isaac   const char *type;
180db4d5e8cSToby Isaac 
181db4d5e8cSToby Isaac   PetscFunctionBegin;
182db4d5e8cSToby Isaac   forest->refct++;
183db4d5e8cSToby Isaac   (*newdm)->data = forest;
1849566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetType((PetscObject)dm, &type));
1859566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, type));
1869566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Forest(*newdm));
1873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
188db4d5e8cSToby Isaac }
189db4d5e8cSToby Isaac 
190d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDestroy_Forest(DM dm)
191d71ae5a4SJacob Faibussowitsch {
192db4d5e8cSToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
193db4d5e8cSToby Isaac 
194db4d5e8cSToby Isaac   PetscFunctionBegin;
1953ba16761SJacob Faibussowitsch   if (--forest->refct > 0) PetscFunctionReturn(PETSC_SUCCESS);
1969566063dSJacob Faibussowitsch   if (forest->destroy) PetscCall((*forest->destroy)(dm));
1979566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&forest->cellSF));
1989566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&forest->preCoarseToFine));
1999566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&forest->coarseToPreFine));
2009566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&forest->adaptLabel));
2019566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest->adaptStrategy));
2029566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&forest->base));
2039566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&forest->adapt));
2049566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest->topology));
2059566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest));
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
207db4d5e8cSToby Isaac }
208db4d5e8cSToby Isaac 
2099be51f97SToby Isaac /*@C
21020f4b53cSBarry Smith   DMForestSetTopology - Set the topology of a `DMFOREST` during the pre-setup phase.  The topology is a string (e.g.
21120f4b53cSBarry Smith   "cube", "shell") and can be interpreted by subtypes of `DMFOREST`) to construct the base DM of a forest during
21220f4b53cSBarry Smith   `DMSetUp()`.
2139be51f97SToby Isaac 
21420f4b53cSBarry Smith   Logically collectiv
2159be51f97SToby Isaac 
21660225df5SJacob Faibussowitsch   Input Parameters:
2179be51f97SToby Isaac + dm       - the forest
2189be51f97SToby Isaac - topology - the topology of the forest
2199be51f97SToby Isaac 
2209be51f97SToby Isaac   Level: intermediate
2219be51f97SToby Isaac 
22220f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetTopology()`, `DMForestSetBaseDM()`
2239be51f97SToby Isaac @*/
224d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology)
225d71ae5a4SJacob Faibussowitsch {
226db4d5e8cSToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
227db4d5e8cSToby Isaac 
228db4d5e8cSToby Isaac   PetscFunctionBegin;
229db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
23028b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the topology after setup");
2319566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest->topology));
2329566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy((const char *)topology, (char **)&forest->topology));
2333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
234db4d5e8cSToby Isaac }
235db4d5e8cSToby Isaac 
2369be51f97SToby Isaac /*@C
23720f4b53cSBarry Smith   DMForestGetTopology - Get a string describing the topology of a `DMFOREST`.
2389be51f97SToby Isaac 
23920f4b53cSBarry Smith   Not Collective
2409be51f97SToby Isaac 
24160225df5SJacob Faibussowitsch   Input Parameter:
2429be51f97SToby Isaac . dm - the forest
2439be51f97SToby Isaac 
24460225df5SJacob Faibussowitsch   Output Parameter:
2459be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell')
2469be51f97SToby Isaac 
2479be51f97SToby Isaac   Level: intermediate
2489be51f97SToby Isaac 
24920f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetTopology()`
2509be51f97SToby Isaac @*/
251d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology)
252d71ae5a4SJacob Faibussowitsch {
253dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
254dd8e54a2SToby Isaac 
255dd8e54a2SToby Isaac   PetscFunctionBegin;
256dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2574f572ea9SToby Isaac   PetscAssertPointer(topology, 2);
258dd8e54a2SToby Isaac   *topology = forest->topology;
2593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
260dd8e54a2SToby Isaac }
261dd8e54a2SToby Isaac 
2629be51f97SToby Isaac /*@
263*a4e35b19SJacob Faibussowitsch   DMForestSetBaseDM - During the pre-setup phase, set the `DM` that defines the base mesh of a
264*a4e35b19SJacob Faibussowitsch   `DMFOREST` forest.
2659be51f97SToby Isaac 
26620f4b53cSBarry Smith   Logically Collective
2679be51f97SToby Isaac 
2689be51f97SToby Isaac   Input Parameters:
2699be51f97SToby Isaac + dm   - the forest
27020f4b53cSBarry Smith - base - the base `DM` of the forest
271765b024eSBarry Smith 
2729be51f97SToby Isaac   Level: intermediate
2739be51f97SToby Isaac 
274*a4e35b19SJacob Faibussowitsch   Notes:
275*a4e35b19SJacob Faibussowitsch   The forest will be hierarchically refined from the base, and all refinements/coarsenings of
276*a4e35b19SJacob Faibussowitsch   the forest will share its base.  In general, two forest must share a base to be comparable,
277*a4e35b19SJacob Faibussowitsch   to do things like construct interpolators.
278*a4e35b19SJacob Faibussowitsch 
27920f4b53cSBarry Smith   Currently the base `DM` must be a `DMPLEX`
28020f4b53cSBarry Smith 
28120f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetBaseDM()`
2829be51f97SToby Isaac @*/
283d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetBaseDM(DM dm, DM base)
284d71ae5a4SJacob Faibussowitsch {
285dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
286dd8e54a2SToby Isaac   PetscInt   dim, dimEmbed;
287dd8e54a2SToby Isaac 
288dd8e54a2SToby Isaac   PetscFunctionBegin;
289dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
29028b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the base after setup");
2919566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)base));
2929566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&forest->base));
293dd8e54a2SToby Isaac   forest->base = base;
294a0452a8eSToby Isaac   if (base) {
2954fb89dddSMatthew G. Knepley     const PetscReal *maxCell, *Lstart, *L;
29628dfcf7cSStefano Zampini 
297a0452a8eSToby Isaac     PetscValidHeaderSpecific(base, DM_CLASSID, 2);
2989566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(base, &dim));
2999566063dSJacob Faibussowitsch     PetscCall(DMSetDimension(dm, dim));
3009566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDim(base, &dimEmbed));
3019566063dSJacob Faibussowitsch     PetscCall(DMSetCoordinateDim(dm, dimEmbed));
3024fb89dddSMatthew G. Knepley     PetscCall(DMGetPeriodicity(base, &maxCell, &Lstart, &L));
3034fb89dddSMatthew G. Knepley     PetscCall(DMSetPeriodicity(dm, maxCell, Lstart, L));
3044fb89dddSMatthew G. Knepley   } else PetscCall(DMSetPeriodicity(dm, NULL, NULL, NULL));
3053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
306dd8e54a2SToby Isaac }
307dd8e54a2SToby Isaac 
3089be51f97SToby Isaac /*@
309*a4e35b19SJacob Faibussowitsch   DMForestGetBaseDM - Get the base `DM` of a `DMFOREST`
3109be51f97SToby Isaac 
31120f4b53cSBarry Smith   Not Collective
3129be51f97SToby Isaac 
3139be51f97SToby Isaac   Input Parameter:
3149be51f97SToby Isaac . dm - the forest
3159be51f97SToby Isaac 
3169be51f97SToby Isaac   Output Parameter:
317*a4e35b19SJacob Faibussowitsch . base - the base `DM` of the forest
318*a4e35b19SJacob Faibussowitsch 
319*a4e35b19SJacob Faibussowitsch   Level: intermediate
3209be51f97SToby Isaac 
321367003a6SStefano Zampini   Notes:
322367003a6SStefano Zampini   After DMSetUp(), the base DM will be redundantly distributed across MPI processes
323367003a6SStefano Zampini 
324*a4e35b19SJacob Faibussowitsch   The forest will be hierarchically refined from the base, and all refinements/coarsenings of
325*a4e35b19SJacob Faibussowitsch   the forest will share its base.  In general, two forest must share a base to be comparable,
326*a4e35b19SJacob Faibussowitsch   to do things like construct interpolators.
3279be51f97SToby Isaac 
32820f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetBaseDM()`
3299be51f97SToby Isaac @*/
330d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetBaseDM(DM dm, DM *base)
331d71ae5a4SJacob Faibussowitsch {
332dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
333dd8e54a2SToby Isaac 
334dd8e54a2SToby Isaac   PetscFunctionBegin;
335dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3364f572ea9SToby Isaac   PetscAssertPointer(base, 2);
337dd8e54a2SToby Isaac   *base = forest->base;
3383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
339dd8e54a2SToby Isaac }
340dd8e54a2SToby Isaac 
341d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *), void *ctx)
342d71ae5a4SJacob Faibussowitsch {
343cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
344cf38a08cSToby Isaac 
345cf38a08cSToby Isaac   PetscFunctionBegin;
346cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
347cf38a08cSToby Isaac   forest->mapcoordinates    = func;
348cf38a08cSToby Isaac   forest->mapcoordinatesctx = ctx;
3493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
350cf38a08cSToby Isaac }
351cf38a08cSToby Isaac 
352d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *), void *ctx)
353d71ae5a4SJacob Faibussowitsch {
354cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
355cf38a08cSToby Isaac 
356cf38a08cSToby Isaac   PetscFunctionBegin;
357cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
358cf38a08cSToby Isaac   if (func) *func = forest->mapcoordinates;
359cf38a08cSToby Isaac   if (ctx) *((void **)ctx) = forest->mapcoordinatesctx;
3603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
361cf38a08cSToby Isaac }
362cf38a08cSToby Isaac 
3639be51f97SToby Isaac /*@
364*a4e35b19SJacob Faibussowitsch   DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the
365*a4e35b19SJacob Faibussowitsch   current forest will be adapted (e.g., the current forest will be
366*a4e35b19SJacob Faibussowitsch   refined/coarsened/repartitioned from it) in `DMSetUp()`.
3679be51f97SToby Isaac 
36820f4b53cSBarry Smith   Logically Collective
3699be51f97SToby Isaac 
370d8d19677SJose E. Roman   Input Parameters:
3719be51f97SToby Isaac + dm    - the new forest, which will be constructed from adapt
3729be51f97SToby Isaac - adapt - the old forest
3739be51f97SToby Isaac 
3749be51f97SToby Isaac   Level: intermediate
3759be51f97SToby Isaac 
37620f4b53cSBarry Smith   Note:
377*a4e35b19SJacob Faibussowitsch   Usually not needed by users directly, `DMForestTemplate()` constructs a new forest to be
378*a4e35b19SJacob Faibussowitsch   adapted from an old forest and calls this routine.
379*a4e35b19SJacob Faibussowitsch 
380*a4e35b19SJacob Faibussowitsch   This can be called after setup with `adapt` = `NULL`, which will clear all internal data
381*a4e35b19SJacob Faibussowitsch   related to the adaptivity forest from `dm`.  This way, repeatedly adapting does not leave
382*a4e35b19SJacob Faibussowitsch   stale `DM` objects in memory.
38320f4b53cSBarry Smith 
38420f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdaptivityForest()`, `DMForestSetAdaptivityPurpose()`
3859be51f97SToby Isaac @*/
386d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityForest(DM dm, DM adapt)
387d71ae5a4SJacob Faibussowitsch {
388dffe73a3SToby Isaac   DM_Forest *forest, *adaptForest, *oldAdaptForest;
389dffe73a3SToby Isaac   DM         oldAdapt;
390456cc5b7SMatthew G. Knepley   PetscBool  isForest;
391dd8e54a2SToby Isaac 
392dd8e54a2SToby Isaac   PetscFunctionBegin;
393dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3941fd50544SStefano Zampini   if (adapt) PetscValidHeaderSpecific(adapt, DM_CLASSID, 2);
3959566063dSJacob Faibussowitsch   PetscCall(DMIsForest(dm, &isForest));
3963ba16761SJacob Faibussowitsch   if (!isForest) PetscFunctionReturn(PETSC_SUCCESS);
3971dca8a05SBarry Smith   PetscCheck(adapt == NULL || !dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adaptation forest after setup");
398ba936b91SToby Isaac   forest = (DM_Forest *)dm->data;
3999566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivityForest(dm, &oldAdapt));
400193eb951SToby Isaac   adaptForest    = (DM_Forest *)(adapt ? adapt->data : NULL);
401193eb951SToby Isaac   oldAdaptForest = (DM_Forest *)(oldAdapt ? oldAdapt->data : NULL);
402dffe73a3SToby Isaac   if (adaptForest != oldAdaptForest) {
4039566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&forest->preCoarseToFine));
4049566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&forest->coarseToPreFine));
4059566063dSJacob Faibussowitsch     if (forest->clearadaptivityforest) PetscCall((*forest->clearadaptivityforest)(dm));
406dffe73a3SToby Isaac   }
40726d9498aSToby Isaac   switch (forest->adaptPurpose) {
408cd3c525cSToby Isaac   case DM_ADAPT_DETERMINE:
4099566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)adapt));
4109566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&(forest->adapt)));
411ba936b91SToby Isaac     forest->adapt = adapt;
41226d9498aSToby Isaac     break;
413d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_REFINE:
414d71ae5a4SJacob Faibussowitsch     PetscCall(DMSetCoarseDM(dm, adapt));
415d71ae5a4SJacob Faibussowitsch     break;
416a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
417d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_COARSEN_LAST:
418d71ae5a4SJacob Faibussowitsch     PetscCall(DMSetFineDM(dm, adapt));
419d71ae5a4SJacob Faibussowitsch     break;
420d71ae5a4SJacob Faibussowitsch   default:
421d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid adaptivity purpose");
42226d9498aSToby Isaac   }
4233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
424dd8e54a2SToby Isaac }
425dd8e54a2SToby Isaac 
4269be51f97SToby Isaac /*@
4279be51f97SToby Isaac   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.
4289be51f97SToby Isaac 
42920f4b53cSBarry Smith   Not Collective
4309be51f97SToby Isaac 
4319be51f97SToby Isaac   Input Parameter:
4329be51f97SToby Isaac . dm - the forest
4339be51f97SToby Isaac 
4349be51f97SToby Isaac   Output Parameter:
43520f4b53cSBarry Smith . adapt - the forest from which `dm` is/was adapted
4369be51f97SToby Isaac 
4379be51f97SToby Isaac   Level: intermediate
4389be51f97SToby Isaac 
43920f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdaptivityForest()`, `DMForestSetAdaptivityPurpose()`
4409be51f97SToby Isaac @*/
441d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
442d71ae5a4SJacob Faibussowitsch {
443ba936b91SToby Isaac   DM_Forest *forest;
444dd8e54a2SToby Isaac 
445dd8e54a2SToby Isaac   PetscFunctionBegin;
446dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
447ba936b91SToby Isaac   forest = (DM_Forest *)dm->data;
44826d9498aSToby Isaac   switch (forest->adaptPurpose) {
449d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_DETERMINE:
450d71ae5a4SJacob Faibussowitsch     *adapt = forest->adapt;
451d71ae5a4SJacob Faibussowitsch     break;
452d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_REFINE:
453d71ae5a4SJacob Faibussowitsch     PetscCall(DMGetCoarseDM(dm, adapt));
454d71ae5a4SJacob Faibussowitsch     break;
455a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
456d71ae5a4SJacob Faibussowitsch   case DM_ADAPT_COARSEN_LAST:
457d71ae5a4SJacob Faibussowitsch     PetscCall(DMGetFineDM(dm, adapt));
458d71ae5a4SJacob Faibussowitsch     break;
459d71ae5a4SJacob Faibussowitsch   default:
460d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid adaptivity purpose");
46126d9498aSToby Isaac   }
4623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46326d9498aSToby Isaac }
46426d9498aSToby Isaac 
4659be51f97SToby Isaac /*@
46620f4b53cSBarry Smith   DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current `DM` is being adapted from its
46720f4b53cSBarry Smith   source (set with `DMForestSetAdaptivityForest()`) for the purpose of refinement (`DM_ADAPT_REFINE`), coarsening
468*a4e35b19SJacob Faibussowitsch   (`DM_ADAPT_COARSEN`), or undefined (`DM_ADAPT_DETERMINE`).
4699be51f97SToby Isaac 
47020f4b53cSBarry Smith   Logically Collective
4719be51f97SToby Isaac 
4729be51f97SToby Isaac   Input Parameters:
4739be51f97SToby Isaac + dm      - the forest
474bf2d5fbbSStefano Zampini - purpose - the adaptivity purpose
4759be51f97SToby Isaac 
4769be51f97SToby Isaac   Level: advanced
4779be51f97SToby Isaac 
478*a4e35b19SJacob Faibussowitsch   Notes:
479*a4e35b19SJacob Faibussowitsch   This only matters for reference counting during `DMDestroy()`. Cyclic references
480*a4e35b19SJacob Faibussowitsch   can be found between `DM`s only if the cyclic reference is due to a fine/coarse relationship
481*a4e35b19SJacob Faibussowitsch   (see `DMSetFineDM()`/`DMSetCoarseDM()`).  If the purpose is not refinement or coarsening, and
482*a4e35b19SJacob Faibussowitsch   the user does not maintain a reference to the post-adaptation forest (i.e., the one created
483*a4e35b19SJacob Faibussowitsch   by `DMForestTemplate()`), this can cause a memory leak.  This method is used by subtypes
484*a4e35b19SJacob Faibussowitsch   of `DMFOREST` when automatically constructing mesh hierarchies.
485*a4e35b19SJacob Faibussowitsch 
48620f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestTemplate()`, `DMForestSetAdaptivityForest()`, `DMForestGetAdaptivityForest()`, `DMAdaptFlag`
4879be51f97SToby Isaac @*/
488d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose)
489d71ae5a4SJacob Faibussowitsch {
49026d9498aSToby Isaac   DM_Forest *forest;
49126d9498aSToby Isaac 
49226d9498aSToby Isaac   PetscFunctionBegin;
49326d9498aSToby Isaac   forest = (DM_Forest *)dm->data;
49428b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adaptation forest after setup");
49526d9498aSToby Isaac   if (purpose != forest->adaptPurpose) {
49626d9498aSToby Isaac     DM adapt;
49726d9498aSToby Isaac 
4989566063dSJacob Faibussowitsch     PetscCall(DMForestGetAdaptivityForest(dm, &adapt));
4999566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)adapt));
5009566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, NULL));
501f885a11aSToby Isaac 
50226d9498aSToby Isaac     forest->adaptPurpose = purpose;
503f885a11aSToby Isaac 
5049566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, adapt));
5059566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&adapt));
50626d9498aSToby Isaac   }
5073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
508dd8e54a2SToby Isaac }
509dd8e54a2SToby Isaac 
51056c0450aSToby Isaac /*@
51120f4b53cSBarry Smith   DMForestGetAdaptivityPurpose - Get whether the current `DM` is being adapted from its source (set with
51220f4b53cSBarry Smith   `DMForestSetAdaptivityForest()`) for the purpose of refinement (`DM_ADAPT_REFINE`), coarsening (`DM_ADAPT_COARSEN`),
51320f4b53cSBarry Smith   coarsening only the last level (`DM_ADAPT_COARSEN_LAST`) or undefined (`DM_ADAPT_DETERMINE`).
51456c0450aSToby Isaac 
51520f4b53cSBarry Smith   Not Collective
51656c0450aSToby Isaac 
51756c0450aSToby Isaac   Input Parameter:
51856c0450aSToby Isaac . dm - the forest
51956c0450aSToby Isaac 
52056c0450aSToby Isaac   Output Parameter:
521bf2d5fbbSStefano Zampini . purpose - the adaptivity purpose
52256c0450aSToby Isaac 
52356c0450aSToby Isaac   Level: advanced
52456c0450aSToby Isaac 
525*a4e35b19SJacob Faibussowitsch   Notes:
526*a4e35b19SJacob Faibussowitsch   This only matters for reference counting: during `DMDestroy()`. Cyclic references
527*a4e35b19SJacob Faibussowitsch   can be found between `DM`s only if the cyclic reference is due to a fine/coarse relationship
528*a4e35b19SJacob Faibussowitsch   (See `DMSetFineDM()`/`DMSetCoarseDM()`).  If the purpose is not refinement or coarsening, and
529*a4e35b19SJacob Faibussowitsch   the user does not maintain a reference to the post-adaptation forest (i.e., the one created
530*a4e35b19SJacob Faibussowitsch   by `DMForestTemplate()`), this can cause a memory leak.  This method is used by subtypes
531*a4e35b19SJacob Faibussowitsch   of `DMFOREST` when automatically constructing mesh hierarchies.
532*a4e35b19SJacob Faibussowitsch 
53320f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestTemplate()`, `DMForestSetAdaptivityForest()`, `DMForestGetAdaptivityForest()`, `DMAdaptFlag`
53456c0450aSToby Isaac @*/
535d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose)
536d71ae5a4SJacob Faibussowitsch {
53756c0450aSToby Isaac   DM_Forest *forest;
53856c0450aSToby Isaac 
53956c0450aSToby Isaac   PetscFunctionBegin;
54056c0450aSToby Isaac   forest   = (DM_Forest *)dm->data;
54156c0450aSToby Isaac   *purpose = forest->adaptPurpose;
5423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
54356c0450aSToby Isaac }
54456c0450aSToby Isaac 
5459be51f97SToby Isaac /*@
5469be51f97SToby Isaac   DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
5479be51f97SToby Isaac   cell adjacency (for the purposes of partitioning and overlap).
5489be51f97SToby Isaac 
54920f4b53cSBarry Smith   Logically Collective
5509be51f97SToby Isaac 
5519be51f97SToby Isaac   Input Parameters:
5529be51f97SToby Isaac + dm     - the forest
5539be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency)
5549be51f97SToby Isaac 
5559be51f97SToby Isaac   Level: intermediate
5569be51f97SToby Isaac 
55720f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`, `DMForestSetPartitionOverlap()`
5589be51f97SToby Isaac @*/
559d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
560d71ae5a4SJacob Faibussowitsch {
561dd8e54a2SToby Isaac   PetscInt   dim;
562dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
563dd8e54a2SToby Isaac 
564dd8e54a2SToby Isaac   PetscFunctionBegin;
565dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
56628b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adjacency dimension after setup");
56763a3b9bcSJacob Faibussowitsch   PetscCheck(adjDim >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "adjacency dim cannot be < 0: %" PetscInt_FMT, adjDim);
5689566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
56963a3b9bcSJacob Faibussowitsch   PetscCheck(adjDim <= dim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "adjacency dim cannot be > %" PetscInt_FMT ": %" PetscInt_FMT, dim, adjDim);
570dd8e54a2SToby Isaac   forest->adjDim = adjDim;
5713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
572dd8e54a2SToby Isaac }
573dd8e54a2SToby Isaac 
5749be51f97SToby Isaac /*@
57520f4b53cSBarry Smith   DMForestSetAdjacencyCodimension - Like `DMForestSetAdjacencyDimension()`, but specified as a co-dimension (so that,
5769be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
5779be51f97SToby Isaac 
57820f4b53cSBarry Smith   Logically Collective
5799be51f97SToby Isaac 
5809be51f97SToby Isaac   Input Parameters:
5819be51f97SToby Isaac + dm       - the forest
58220f4b53cSBarry Smith - adjCodim - default is the dimension of the forest (see `DMGetDimension()`), since this is the codimension of vertices
5839be51f97SToby Isaac 
5849be51f97SToby Isaac   Level: intermediate
5859be51f97SToby Isaac 
58620f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdjacencyCodimension()`, `DMForestSetAdjacencyDimension()`
5879be51f97SToby Isaac @*/
588d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
589d71ae5a4SJacob Faibussowitsch {
590dd8e54a2SToby Isaac   PetscInt dim;
591dd8e54a2SToby Isaac 
592dd8e54a2SToby Isaac   PetscFunctionBegin;
593dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5949566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
5959566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdjacencyDimension(dm, dim - adjCodim));
5963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
597dd8e54a2SToby Isaac }
598dd8e54a2SToby Isaac 
5999be51f97SToby Isaac /*@
6009be51f97SToby Isaac   DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
6019be51f97SToby Isaac   purposes of partitioning and overlap).
6029be51f97SToby Isaac 
60320f4b53cSBarry Smith   Not Collective
6049be51f97SToby Isaac 
6059be51f97SToby Isaac   Input Parameter:
6069be51f97SToby Isaac . dm - the forest
6079be51f97SToby Isaac 
6089be51f97SToby Isaac   Output Parameter:
6099be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency)
6109be51f97SToby Isaac 
6119be51f97SToby Isaac   Level: intermediate
6129be51f97SToby Isaac 
61320f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdjacencyDimension()`, `DMForestGetAdjacencyCodimension()`, `DMForestSetPartitionOverlap()`
6149be51f97SToby Isaac @*/
615d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
616d71ae5a4SJacob Faibussowitsch {
617dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
618dd8e54a2SToby Isaac 
619dd8e54a2SToby Isaac   PetscFunctionBegin;
620dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6214f572ea9SToby Isaac   PetscAssertPointer(adjDim, 2);
622dd8e54a2SToby Isaac   *adjDim = forest->adjDim;
6233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
624dd8e54a2SToby Isaac }
625dd8e54a2SToby Isaac 
6269be51f97SToby Isaac /*@
62720f4b53cSBarry Smith   DMForestGetAdjacencyCodimension - Like `DMForestGetAdjacencyDimension()`, but specified as a co-dimension (so that,
6289be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
6299be51f97SToby Isaac 
63020f4b53cSBarry Smith   Not Collective
6319be51f97SToby Isaac 
6329be51f97SToby Isaac   Input Parameter:
6339be51f97SToby Isaac . dm - the forest
6349be51f97SToby Isaac 
6359be51f97SToby Isaac   Output Parameter:
63620f4b53cSBarry Smith . adjCodim - default isthe dimension of the forest (see `DMGetDimension()`), since this is the codimension of vertices
6379be51f97SToby Isaac 
6389be51f97SToby Isaac   Level: intermediate
6399be51f97SToby Isaac 
64020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdjacencyCodimension()`, `DMForestGetAdjacencyDimension()`
6419be51f97SToby Isaac @*/
642d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
643d71ae5a4SJacob Faibussowitsch {
644dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
645dd8e54a2SToby Isaac   PetscInt   dim;
646dd8e54a2SToby Isaac 
647dd8e54a2SToby Isaac   PetscFunctionBegin;
648dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
6494f572ea9SToby Isaac   PetscAssertPointer(adjCodim, 2);
6509566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
651dd8e54a2SToby Isaac   *adjCodim = dim - forest->adjDim;
6523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
653dd8e54a2SToby Isaac }
654dd8e54a2SToby Isaac 
6559be51f97SToby Isaac /*@
6569be51f97SToby Isaac   DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel
6579be51f97SToby Isaac   partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding
6589be51f97SToby Isaac   adjacent cells
6599be51f97SToby Isaac 
66020f4b53cSBarry Smith   Logically Collective
6619be51f97SToby Isaac 
6629be51f97SToby Isaac   Input Parameters:
6639be51f97SToby Isaac + dm      - the forest
6649be51f97SToby Isaac - overlap - default 0
6659be51f97SToby Isaac 
6669be51f97SToby Isaac   Level: intermediate
6679be51f97SToby Isaac 
66820f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetPartitionOverlap()`, `DMForestSetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`
6699be51f97SToby Isaac @*/
670d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap)
671d71ae5a4SJacob Faibussowitsch {
672dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
673dd8e54a2SToby Isaac 
674dd8e54a2SToby Isaac   PetscFunctionBegin;
675dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
67628b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the overlap after setup");
67763a3b9bcSJacob Faibussowitsch   PetscCheck(overlap >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "overlap cannot be < 0: %" PetscInt_FMT, overlap);
678dd8e54a2SToby Isaac   forest->overlap = overlap;
6793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
680dd8e54a2SToby Isaac }
681dd8e54a2SToby Isaac 
6829be51f97SToby Isaac /*@
6839be51f97SToby Isaac   DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values
6849be51f97SToby Isaac   > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells
6859be51f97SToby Isaac 
68620f4b53cSBarry Smith   Not Collective
6879be51f97SToby Isaac 
6889be51f97SToby Isaac   Input Parameter:
6899be51f97SToby Isaac . dm - the forest
6909be51f97SToby Isaac 
6919be51f97SToby Isaac   Output Parameter:
6929be51f97SToby Isaac . overlap - default 0
6939be51f97SToby Isaac 
6949be51f97SToby Isaac   Level: intermediate
6959be51f97SToby Isaac 
69642747ad1SJacob Faibussowitsch .seealso: `DM`, `DMFOREST`, `DMForestSetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`
6979be51f97SToby Isaac @*/
698d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap)
699d71ae5a4SJacob Faibussowitsch {
700dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
701dd8e54a2SToby Isaac 
702dd8e54a2SToby Isaac   PetscFunctionBegin;
703dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7044f572ea9SToby Isaac   PetscAssertPointer(overlap, 2);
705dd8e54a2SToby Isaac   *overlap = forest->overlap;
7063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
707dd8e54a2SToby Isaac }
708dd8e54a2SToby Isaac 
7099be51f97SToby Isaac /*@
7109be51f97SToby Isaac   DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base
71120f4b53cSBarry Smith   `DM`, see `DMForestGetBaseDM()`) allowed in the forest.  If the forest is being created by coarsening a previous forest
71220f4b53cSBarry Smith   (see `DMForestGetAdaptivityForest()`) this limits the amount of coarsening.
7139be51f97SToby Isaac 
71420f4b53cSBarry Smith   Logically Collective
7159be51f97SToby Isaac 
7169be51f97SToby Isaac   Input Parameters:
7179be51f97SToby Isaac + dm            - the forest
71820f4b53cSBarry Smith - minRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`)
7199be51f97SToby Isaac 
7209be51f97SToby Isaac   Level: intermediate
7219be51f97SToby Isaac 
72220f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestSetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()`
7239be51f97SToby Isaac @*/
724d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement)
725d71ae5a4SJacob Faibussowitsch {
726dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
727dd8e54a2SToby Isaac 
728dd8e54a2SToby Isaac   PetscFunctionBegin;
729dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
73028b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the minimum refinement after setup");
731dd8e54a2SToby Isaac   forest->minRefinement = minRefinement;
7323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
733dd8e54a2SToby Isaac }
734dd8e54a2SToby Isaac 
7359be51f97SToby Isaac /*@
73620f4b53cSBarry Smith   DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base `DM`, see
73720f4b53cSBarry Smith   `DMForestGetBaseDM()`) allowed in the forest.  If the forest is being created by coarsening a previous forest (see
73820f4b53cSBarry Smith   `DMForestGetAdaptivityForest()`), this limits the amount of coarsening.
7399be51f97SToby Isaac 
74020f4b53cSBarry Smith   Not Collective
7419be51f97SToby Isaac 
7429be51f97SToby Isaac   Input Parameter:
7439be51f97SToby Isaac . dm - the forest
7449be51f97SToby Isaac 
7459be51f97SToby Isaac   Output Parameter:
74620f4b53cSBarry Smith . minRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`)
7479be51f97SToby Isaac 
7489be51f97SToby Isaac   Level: intermediate
7499be51f97SToby Isaac 
75020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetMinimumRefinement()`, `DMForestGetMaximumRefinement()`, `DMForestGetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()`
7519be51f97SToby Isaac @*/
752d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement)
753d71ae5a4SJacob Faibussowitsch {
754dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
755dd8e54a2SToby Isaac 
756dd8e54a2SToby Isaac   PetscFunctionBegin;
757dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7584f572ea9SToby Isaac   PetscAssertPointer(minRefinement, 2);
759dd8e54a2SToby Isaac   *minRefinement = forest->minRefinement;
7603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
761dd8e54a2SToby Isaac }
762dd8e54a2SToby Isaac 
7639be51f97SToby Isaac /*@
7649be51f97SToby Isaac   DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base
76520f4b53cSBarry Smith   `DM`, see `DMForestGetBaseDM()`) allowed in the forest.
7669be51f97SToby Isaac 
76720f4b53cSBarry Smith   Logically Collective
7689be51f97SToby Isaac 
7699be51f97SToby Isaac   Input Parameters:
7709be51f97SToby Isaac + dm             - the forest
77160225df5SJacob Faibussowitsch - initRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`)
7729be51f97SToby Isaac 
7739be51f97SToby Isaac   Level: intermediate
7749be51f97SToby Isaac 
77520f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestGetBaseDM()`
7769be51f97SToby Isaac @*/
777d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement)
778d71ae5a4SJacob Faibussowitsch {
77956ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
78056ba9f64SToby Isaac 
78156ba9f64SToby Isaac   PetscFunctionBegin;
78256ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
78328b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the initial refinement after setup");
78456ba9f64SToby Isaac   forest->initRefinement = initRefinement;
7853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
78656ba9f64SToby Isaac }
78756ba9f64SToby Isaac 
7889be51f97SToby Isaac /*@
78920f4b53cSBarry Smith   DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base `DM`, see
79020f4b53cSBarry Smith   `DMForestGetBaseDM()`) allowed in the forest.
7919be51f97SToby Isaac 
79220f4b53cSBarry Smith   Not Collective
7939be51f97SToby Isaac 
7949be51f97SToby Isaac   Input Parameter:
7959be51f97SToby Isaac . dm - the forest
7969be51f97SToby Isaac 
79701d2d390SJose E. Roman   Output Parameter:
79820f4b53cSBarry Smith . initRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`)
7999be51f97SToby Isaac 
8009be51f97SToby Isaac   Level: intermediate
8019be51f97SToby Isaac 
80220f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestGetBaseDM()`
8039be51f97SToby Isaac @*/
804d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement)
805d71ae5a4SJacob Faibussowitsch {
80656ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
80756ba9f64SToby Isaac 
80856ba9f64SToby Isaac   PetscFunctionBegin;
80956ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8104f572ea9SToby Isaac   PetscAssertPointer(initRefinement, 2);
81156ba9f64SToby Isaac   *initRefinement = forest->initRefinement;
8123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
81356ba9f64SToby Isaac }
81456ba9f64SToby Isaac 
8159be51f97SToby Isaac /*@
8169be51f97SToby Isaac   DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base
81720f4b53cSBarry Smith   `DM`, see `DMForestGetBaseDM()`) allowed in the forest.  If the forest is being created by refining a previous forest
81820f4b53cSBarry Smith   (see `DMForestGetAdaptivityForest()`), this limits the amount of refinement.
8199be51f97SToby Isaac 
82020f4b53cSBarry Smith   Logically Collective
8219be51f97SToby Isaac 
8229be51f97SToby Isaac   Input Parameters:
8239be51f97SToby Isaac + dm            - the forest
82420f4b53cSBarry Smith - maxRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`)
8259be51f97SToby Isaac 
8269be51f97SToby Isaac   Level: intermediate
8279be51f97SToby Isaac 
82842747ad1SJacob Faibussowitsch .seealso: `DM`, `DMFOREST`, `DMForestGetMinimumRefinement()`, `DMForestSetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityDM()`
8299be51f97SToby Isaac @*/
830d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement)
831d71ae5a4SJacob Faibussowitsch {
832dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
833dd8e54a2SToby Isaac 
834dd8e54a2SToby Isaac   PetscFunctionBegin;
835dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
83628b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the maximum refinement after setup");
837c7eeac06SToby Isaac   forest->maxRefinement = maxRefinement;
8383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
839dd8e54a2SToby Isaac }
840dd8e54a2SToby Isaac 
8419be51f97SToby Isaac /*@
84220f4b53cSBarry Smith   DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base `DM`, see
84320f4b53cSBarry Smith   `DMForestGetBaseDM()`) allowed in the forest.  If the forest is being created by refining a previous forest (see
844f826b5fcSPierre Jolivet   `DMForestGetAdaptivityForest()`), this limits the amount of refinement.
8459be51f97SToby Isaac 
84620f4b53cSBarry Smith   Not Collective
8479be51f97SToby Isaac 
8489be51f97SToby Isaac   Input Parameter:
8499be51f97SToby Isaac . dm - the forest
8509be51f97SToby Isaac 
8519be51f97SToby Isaac   Output Parameter:
85220f4b53cSBarry Smith . maxRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`)
8539be51f97SToby Isaac 
8549be51f97SToby Isaac   Level: intermediate
8559be51f97SToby Isaac 
85620f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetMaximumRefinement()`, `DMForestGetMinimumRefinement()`, `DMForestGetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()`
8579be51f97SToby Isaac @*/
858d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement)
859d71ae5a4SJacob Faibussowitsch {
860dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
861dd8e54a2SToby Isaac 
862dd8e54a2SToby Isaac   PetscFunctionBegin;
863dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8644f572ea9SToby Isaac   PetscAssertPointer(maxRefinement, 2);
865c7eeac06SToby Isaac   *maxRefinement = forest->maxRefinement;
8663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
867dd8e54a2SToby Isaac }
868c7eeac06SToby Isaac 
8699be51f97SToby Isaac /*@C
8709be51f97SToby Isaac   DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes.
8719be51f97SToby Isaac 
87220f4b53cSBarry Smith   Logically Collective
8739be51f97SToby Isaac 
8749be51f97SToby Isaac   Input Parameters:
8759be51f97SToby Isaac + dm            - the forest
87620f4b53cSBarry Smith - adaptStrategy - default `DMFORESTADAPTALL`
8779be51f97SToby Isaac 
8789be51f97SToby Isaac   Level: advanced
8799be51f97SToby Isaac 
88020f4b53cSBarry Smith   Notes:
88120f4b53cSBarry Smith   Subtypes of `DMFOREST` may define their own strategies.  Two default strategies are `DMFORESTADAPTALL`, which indicates that all processes must agree
88220f4b53cSBarry Smith   for a refinement/coarsening flag to be valid, and `DMFORESTADAPTANY`, which indicates that only one process needs to
88320f4b53cSBarry Smith   specify refinement/coarsening.
88420f4b53cSBarry Smith 
88520f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdaptivityStrategy()`, `DMFORESTADAPTALL`, `DMFORESTADAPTANY`
8869be51f97SToby Isaac @*/
887d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy)
888d71ae5a4SJacob Faibussowitsch {
889c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
890c7eeac06SToby Isaac 
891c7eeac06SToby Isaac   PetscFunctionBegin;
892c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8939566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest->adaptStrategy));
8949566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy((const char *)adaptStrategy, (char **)&forest->adaptStrategy));
8953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
896c7eeac06SToby Isaac }
897c7eeac06SToby Isaac 
8989be51f97SToby Isaac /*@C
89960225df5SJacob Faibussowitsch   DMForestGetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes.
9009be51f97SToby Isaac 
90120f4b53cSBarry Smith   Not Collective
9029be51f97SToby Isaac 
9039be51f97SToby Isaac   Input Parameter:
9049be51f97SToby Isaac . dm - the forest
9059be51f97SToby Isaac 
9069be51f97SToby Isaac   Output Parameter:
90720f4b53cSBarry Smith . adaptStrategy - the adaptivity strategy (default `DMFORESTADAPTALL`)
9089be51f97SToby Isaac 
9099be51f97SToby Isaac   Level: advanced
9109be51f97SToby Isaac 
91120f4b53cSBarry Smith   Note:
91220f4b53cSBarry Smith   Subtypes
91320f4b53cSBarry Smith   of `DMFOREST` may define their own strategies.  Two default strategies are `DMFORESTADAPTALL`, which indicates that all
91420f4b53cSBarry Smith   processes must agree for a refinement/coarsening flag to be valid, and `DMFORESTADAPTANY`, which indicates that only
91520f4b53cSBarry Smith   one process needs to specify refinement/coarsening.
91620f4b53cSBarry Smith 
91720f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMFORESTADAPTALL`, `DMFORESTADAPTANY`, `DMForestSetAdaptivityStrategy()`
9189be51f97SToby Isaac @*/
919d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
920d71ae5a4SJacob Faibussowitsch {
921c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
922c7eeac06SToby Isaac 
923c7eeac06SToby Isaac   PetscFunctionBegin;
924c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9254f572ea9SToby Isaac   PetscAssertPointer(adaptStrategy, 2);
926c7eeac06SToby Isaac   *adaptStrategy = forest->adaptStrategy;
9273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
928c7eeac06SToby Isaac }
929c7eeac06SToby Isaac 
9302a133e43SToby Isaac /*@
9312a133e43SToby Isaac   DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning,
93220f4b53cSBarry Smith   etc.) was successful.
9332a133e43SToby Isaac 
93420f4b53cSBarry Smith   Collective
9352a133e43SToby Isaac 
9362a133e43SToby Isaac   Input Parameter:
9372a133e43SToby Isaac . dm - the post-adaptation forest
9382a133e43SToby Isaac 
9392a133e43SToby Isaac   Output Parameter:
94020f4b53cSBarry Smith . success - `PETSC_TRUE` if the post-adaptation forest is different from the pre-adaptation forest.
9412a133e43SToby Isaac 
9422a133e43SToby Isaac   Level: intermediate
9432a133e43SToby Isaac 
94420f4b53cSBarry Smith   Notes:
94520f4b53cSBarry Smith   `PETSC_FALSE` indicates that the post-adaptation forest is the same as the pre-adpatation
94620f4b53cSBarry Smith   forest.  A requested adaptation may have been unsuccessful if, for example, the requested refinement would have
94720f4b53cSBarry Smith   exceeded the maximum refinement level.
94820f4b53cSBarry Smith 
94920f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`
9502a133e43SToby Isaac @*/
951d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success)
952d71ae5a4SJacob Faibussowitsch {
9532a133e43SToby Isaac   DM_Forest *forest;
9542a133e43SToby Isaac 
9552a133e43SToby Isaac   PetscFunctionBegin;
9562a133e43SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
95728b400f6SJacob Faibussowitsch   PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() has not been called yet.");
9582a133e43SToby Isaac   forest = (DM_Forest *)dm->data;
9599566063dSJacob Faibussowitsch   PetscCall((forest->getadaptivitysuccess)(dm, success));
9603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9612a133e43SToby Isaac }
9622a133e43SToby Isaac 
963bf9b5d84SToby Isaac /*@
96420f4b53cSBarry Smith   DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer `PetscSF`s should be computed
96520f4b53cSBarry Smith   relating the cells of the pre-adaptation forest to the post-adaptiation forest.
966bf9b5d84SToby Isaac 
96720f4b53cSBarry Smith   Logically Collective
968bf9b5d84SToby Isaac 
969bf9b5d84SToby Isaac   Input Parameters:
970bf9b5d84SToby Isaac + dm        - the post-adaptation forest
97120f4b53cSBarry Smith - computeSF - default `PETSC_TRUE`
972bf9b5d84SToby Isaac 
973bf9b5d84SToby Isaac   Level: advanced
974bf9b5d84SToby Isaac 
97520f4b53cSBarry Smith   Note:
97620f4b53cSBarry Smith   After `DMSetUp()` is called, the transfer `PetscSF`s can be accessed with `DMForestGetAdaptivitySF()`.
97720f4b53cSBarry Smith 
97820f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetComputeAdaptivitySF()`, `DMForestGetAdaptivitySF()`
979bf9b5d84SToby Isaac @*/
980d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
981d71ae5a4SJacob Faibussowitsch {
982bf9b5d84SToby Isaac   DM_Forest *forest;
983bf9b5d84SToby Isaac 
984bf9b5d84SToby Isaac   PetscFunctionBegin;
985bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
98628b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot compute adaptivity PetscSFs after setup is called");
987bf9b5d84SToby Isaac   forest                 = (DM_Forest *)dm->data;
988bf9b5d84SToby Isaac   forest->computeAdaptSF = computeSF;
9893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
990bf9b5d84SToby Isaac }
991bf9b5d84SToby Isaac 
992d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time)
993d71ae5a4SJacob Faibussowitsch {
99480b27e07SToby Isaac   DM_Forest *forest;
99580b27e07SToby Isaac 
99680b27e07SToby Isaac   PetscFunctionBegin;
99780b27e07SToby Isaac   PetscValidHeaderSpecific(dmIn, DM_CLASSID, 1);
99880b27e07SToby Isaac   PetscValidHeaderSpecific(vecIn, VEC_CLASSID, 2);
99980b27e07SToby Isaac   PetscValidHeaderSpecific(dmOut, DM_CLASSID, 3);
100080b27e07SToby Isaac   PetscValidHeaderSpecific(vecOut, VEC_CLASSID, 4);
100180b27e07SToby Isaac   forest = (DM_Forest *)dmIn->data;
100228b400f6SJacob Faibussowitsch   PetscCheck(forest->transfervec, PetscObjectComm((PetscObject)dmIn), PETSC_ERR_SUP, "DMForestTransferVec() not implemented");
10039566063dSJacob Faibussowitsch   PetscCall((forest->transfervec)(dmIn, vecIn, dmOut, vecOut, useBCs, time));
10043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
100580b27e07SToby Isaac }
100680b27e07SToby Isaac 
1007d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestTransferVecFromBase(DM dm, Vec vecIn, Vec vecOut)
1008d71ae5a4SJacob Faibussowitsch {
1009ac34a06fSStefano Zampini   DM_Forest *forest;
1010ac34a06fSStefano Zampini 
1011ac34a06fSStefano Zampini   PetscFunctionBegin;
1012ac34a06fSStefano Zampini   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1013ac34a06fSStefano Zampini   PetscValidHeaderSpecific(vecIn, VEC_CLASSID, 2);
1014ac34a06fSStefano Zampini   PetscValidHeaderSpecific(vecOut, VEC_CLASSID, 3);
1015ac34a06fSStefano Zampini   forest = (DM_Forest *)dm->data;
101628b400f6SJacob Faibussowitsch   PetscCheck(forest->transfervecfrombase, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMForestTransferVecFromBase() not implemented");
10179566063dSJacob Faibussowitsch   PetscCall((forest->transfervecfrombase)(dm, vecIn, vecOut));
10183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1019ac34a06fSStefano Zampini }
1020ac34a06fSStefano Zampini 
1021bf9b5d84SToby Isaac /*@
102220f4b53cSBarry Smith   DMForestGetComputeAdaptivitySF - Get whether transfer `PetscSF`s should be computed relating the cells of the
102320f4b53cSBarry Smith   pre-adaptation forest to the post-adaptiation forest.  After `DMSetUp()` is called, these transfer PetscSFs can be
102420f4b53cSBarry Smith   accessed with `DMForestGetAdaptivitySF()`.
1025bf9b5d84SToby Isaac 
102620f4b53cSBarry Smith   Not Collective
1027bf9b5d84SToby Isaac 
1028bf9b5d84SToby Isaac   Input Parameter:
1029bf9b5d84SToby Isaac . dm - the post-adaptation forest
1030bf9b5d84SToby Isaac 
1031bf9b5d84SToby Isaac   Output Parameter:
103220f4b53cSBarry Smith . computeSF - default `PETSC_TRUE`
1033bf9b5d84SToby Isaac 
1034bf9b5d84SToby Isaac   Level: advanced
1035bf9b5d84SToby Isaac 
103620f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetComputeAdaptivitySF()`, `DMForestGetAdaptivitySF()`
1037bf9b5d84SToby Isaac @*/
1038d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1039d71ae5a4SJacob Faibussowitsch {
1040bf9b5d84SToby Isaac   DM_Forest *forest;
1041bf9b5d84SToby Isaac 
1042bf9b5d84SToby Isaac   PetscFunctionBegin;
1043bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1044bf9b5d84SToby Isaac   forest     = (DM_Forest *)dm->data;
1045bf9b5d84SToby Isaac   *computeSF = forest->computeAdaptSF;
10463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1047bf9b5d84SToby Isaac }
1048bf9b5d84SToby Isaac 
1049bf9b5d84SToby Isaac /*@
1050*a4e35b19SJacob Faibussowitsch   DMForestGetAdaptivitySF - Get `PetscSF`s that relate the pre-adaptation forest to the
1051*a4e35b19SJacob Faibussowitsch   post-adaptation forest.
1052bf9b5d84SToby Isaac 
105320f4b53cSBarry Smith   Not Collective
1054bf9b5d84SToby Isaac 
1055bf9b5d84SToby Isaac   Input Parameter:
105620f4b53cSBarry Smith . dm - the post-adaptation forest
1057bf9b5d84SToby Isaac 
105820f4b53cSBarry Smith   Output Parameters:
105920f4b53cSBarry Smith + preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
106020f4b53cSBarry Smith - coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-
1061bf9b5d84SToby Isaac 
1062bf9b5d84SToby Isaac   Level: advanced
1063bf9b5d84SToby Isaac 
1064*a4e35b19SJacob Faibussowitsch   Notes:
1065*a4e35b19SJacob Faibussowitsch   Adaptation can be any combination of refinement, coarsening, repartition, and change of
1066*a4e35b19SJacob Faibussowitsch   overlap, so there may be some cells of the pre-adaptation that are parents of post-adaptation
1067*a4e35b19SJacob Faibussowitsch   cells, and vice versa.  Therefore there are two `PetscSF`s: one that relates pre-adaptation
1068*a4e35b19SJacob Faibussowitsch   coarse cells to post-adaptation fine cells, and one that relates pre-adaptation fine cells to
1069*a4e35b19SJacob Faibussowitsch   post-adaptation coarse cells.
1070*a4e35b19SJacob Faibussowitsch 
107120f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetComputeAdaptivitySF()`, `DMForestSetComputeAdaptivitySF()`
1072bf9b5d84SToby Isaac @*/
1073d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1074d71ae5a4SJacob Faibussowitsch {
1075bf9b5d84SToby Isaac   DM_Forest *forest;
1076bf9b5d84SToby Isaac 
1077bf9b5d84SToby Isaac   PetscFunctionBegin;
1078bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10799566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
1080bf9b5d84SToby Isaac   forest = (DM_Forest *)dm->data;
1081f885a11aSToby Isaac   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1082f885a11aSToby Isaac   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
10833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1084bf9b5d84SToby Isaac }
1085bf9b5d84SToby Isaac 
10869be51f97SToby Isaac /*@
1087*a4e35b19SJacob Faibussowitsch   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the
1088*a4e35b19SJacob Faibussowitsch   mesh, e.g. give 2 to indicate that the diameter of neighboring cells should differ by at most
1089*a4e35b19SJacob Faibussowitsch   a factor of 2.
10909be51f97SToby Isaac 
109120f4b53cSBarry Smith   Logically Collective
10929be51f97SToby Isaac 
10939be51f97SToby Isaac   Input Parameters:
10949be51f97SToby Isaac + dm    - the forest
10959be51f97SToby Isaac - grade - the grading factor
10969be51f97SToby Isaac 
10979be51f97SToby Isaac   Level: advanced
10989be51f97SToby Isaac 
1099*a4e35b19SJacob Faibussowitsch   Note:
1100*a4e35b19SJacob Faibussowitsch   Subtypes of `DMFOREST` may only support one particular choice of grading factor.
1101*a4e35b19SJacob Faibussowitsch 
110220f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetGradeFactor()`
11039be51f97SToby Isaac @*/
1104d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1105d71ae5a4SJacob Faibussowitsch {
1106c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1107c7eeac06SToby Isaac 
1108c7eeac06SToby Isaac   PetscFunctionBegin;
1109c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
111028b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the grade factor after setup");
1111c7eeac06SToby Isaac   forest->gradeFactor = grade;
11123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1113c7eeac06SToby Isaac }
1114c7eeac06SToby Isaac 
11159be51f97SToby Isaac /*@
11169be51f97SToby Isaac   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
111720f4b53cSBarry Smith   neighboring cells should differ by at most a factor of 2.  Subtypes of `DMFOREST` may only support one particular
11189be51f97SToby Isaac   choice of grading factor.
11199be51f97SToby Isaac 
112020f4b53cSBarry Smith   Not Collective
11219be51f97SToby Isaac 
11229be51f97SToby Isaac   Input Parameter:
11239be51f97SToby Isaac . dm - the forest
11249be51f97SToby Isaac 
11259be51f97SToby Isaac   Output Parameter:
11269be51f97SToby Isaac . grade - the grading factor
11279be51f97SToby Isaac 
11289be51f97SToby Isaac   Level: advanced
11299be51f97SToby Isaac 
113020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetGradeFactor()`
11319be51f97SToby Isaac @*/
1132d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1133d71ae5a4SJacob Faibussowitsch {
1134c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1135c7eeac06SToby Isaac 
1136c7eeac06SToby Isaac   PetscFunctionBegin;
1137c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11384f572ea9SToby Isaac   PetscAssertPointer(grade, 2);
1139c7eeac06SToby Isaac   *grade = forest->gradeFactor;
11403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1141c7eeac06SToby Isaac }
1142c7eeac06SToby Isaac 
11439be51f97SToby Isaac /*@
11449be51f97SToby Isaac   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
1145*a4e35b19SJacob Faibussowitsch   the cell weight (see `DMForestSetCellWeights()`) when calculating partitions.
11469be51f97SToby Isaac 
114720f4b53cSBarry Smith   Logically Collective
11489be51f97SToby Isaac 
11499be51f97SToby Isaac   Input Parameters:
11509be51f97SToby Isaac + dm            - the forest
115160225df5SJacob Faibussowitsch - weightsFactor - default 1.
11529be51f97SToby Isaac 
11539be51f97SToby Isaac   Level: advanced
11549be51f97SToby Isaac 
1155*a4e35b19SJacob Faibussowitsch   Note:
1156*a4e35b19SJacob Faibussowitsch   The final weight of a cell will be (cellWeight) * (weightFactor^refinementLevel).  A factor
1157*a4e35b19SJacob Faibussowitsch   of 1 indicates that the weight of a cell does not depend on its level; a factor of 2, for
1158*a4e35b19SJacob Faibussowitsch   example, might be appropriate for sub-cycling time-stepping methods, when the computation
1159*a4e35b19SJacob Faibussowitsch   associated with a cell is multiplied by a factor of 2 for each additional level of
1160*a4e35b19SJacob Faibussowitsch   refinement.
1161*a4e35b19SJacob Faibussowitsch 
116220f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetCellWeightFactor()`, `DMForestSetCellWeights()`
11639be51f97SToby Isaac @*/
1164d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
1165d71ae5a4SJacob Faibussowitsch {
1166c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1167c7eeac06SToby Isaac 
1168c7eeac06SToby Isaac   PetscFunctionBegin;
1169c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
117028b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the weights factor after setup");
1171c7eeac06SToby Isaac   forest->weightsFactor = weightsFactor;
11723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1173c7eeac06SToby Isaac }
1174c7eeac06SToby Isaac 
11759be51f97SToby Isaac /*@
11769be51f97SToby Isaac   DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
1177*a4e35b19SJacob Faibussowitsch   `DMForestSetCellWeights()`) when calculating partitions.
11789be51f97SToby Isaac 
117920f4b53cSBarry Smith   Not Collective
11809be51f97SToby Isaac 
11819be51f97SToby Isaac   Input Parameter:
11829be51f97SToby Isaac . dm - the forest
11839be51f97SToby Isaac 
11849be51f97SToby Isaac   Output Parameter:
118560225df5SJacob Faibussowitsch . weightsFactor - default 1.
11869be51f97SToby Isaac 
11879be51f97SToby Isaac   Level: advanced
11889be51f97SToby Isaac 
1189*a4e35b19SJacob Faibussowitsch   Note:
1190*a4e35b19SJacob Faibussowitsch   The final weight of a cell will be (cellWeight) * (weightFactor^refinementLevel).  A factor
1191*a4e35b19SJacob Faibussowitsch   of 1 indicates that the weight of a cell does not depend on its level; a factor of 2, for
1192*a4e35b19SJacob Faibussowitsch   example, might be appropriate for sub-cycling time-stepping methods, when the computation
1193*a4e35b19SJacob Faibussowitsch   associated with a cell is multiplied by a factor of 2 for each additional level of
1194*a4e35b19SJacob Faibussowitsch   refinement.
1195*a4e35b19SJacob Faibussowitsch 
119620f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetCellWeightFactor()`, `DMForestSetCellWeights()`
11979be51f97SToby Isaac @*/
1198d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1199d71ae5a4SJacob Faibussowitsch {
1200c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1201c7eeac06SToby Isaac 
1202c7eeac06SToby Isaac   PetscFunctionBegin;
1203c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12044f572ea9SToby Isaac   PetscAssertPointer(weightsFactor, 2);
1205c7eeac06SToby Isaac   *weightsFactor = forest->weightsFactor;
12063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1207c7eeac06SToby Isaac }
1208c7eeac06SToby Isaac 
12099be51f97SToby Isaac /*@
12109be51f97SToby Isaac   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
12119be51f97SToby Isaac 
121220f4b53cSBarry Smith   Not Collective
12139be51f97SToby Isaac 
12149be51f97SToby Isaac   Input Parameter:
12159be51f97SToby Isaac . dm - the forest
12169be51f97SToby Isaac 
12179be51f97SToby Isaac   Output Parameters:
12189be51f97SToby Isaac + cStart - the first cell on this process
12199be51f97SToby Isaac - cEnd   - one after the final cell on this process
12209be51f97SToby Isaac 
12211a244344SSatish Balay   Level: intermediate
12229be51f97SToby Isaac 
122320f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetCellSF()`
12249be51f97SToby Isaac @*/
1225d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1226d71ae5a4SJacob Faibussowitsch {
1227c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1228c7eeac06SToby Isaac 
1229c7eeac06SToby Isaac   PetscFunctionBegin;
1230c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12314f572ea9SToby Isaac   PetscAssertPointer(cStart, 2);
12324f572ea9SToby Isaac   PetscAssertPointer(cEnd, 3);
123348a46eb9SPierre Jolivet   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) PetscCall(forest->createcellchart(dm, &forest->cStart, &forest->cEnd));
1234c7eeac06SToby Isaac   *cStart = forest->cStart;
1235c7eeac06SToby Isaac   *cEnd   = forest->cEnd;
12363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1237c7eeac06SToby Isaac }
1238c7eeac06SToby Isaac 
12399be51f97SToby Isaac /*@
124020f4b53cSBarry Smith   DMForestGetCellSF - After the setup phase, get the `PetscSF` for overlapping cells between processes
12419be51f97SToby Isaac 
124220f4b53cSBarry Smith   Not Collective
12439be51f97SToby Isaac 
12449be51f97SToby Isaac   Input Parameter:
12459be51f97SToby Isaac . dm - the forest
12469be51f97SToby Isaac 
12479be51f97SToby Isaac   Output Parameter:
124820f4b53cSBarry Smith . cellSF - the `PetscSF`
12499be51f97SToby Isaac 
12501a244344SSatish Balay   Level: intermediate
12519be51f97SToby Isaac 
125220f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetCellChart()`
12539be51f97SToby Isaac @*/
1254d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1255d71ae5a4SJacob Faibussowitsch {
1256c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1257c7eeac06SToby Isaac 
1258c7eeac06SToby Isaac   PetscFunctionBegin;
1259c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12604f572ea9SToby Isaac   PetscAssertPointer(cellSF, 2);
126148a46eb9SPierre Jolivet   if ((!forest->cellSF) && forest->createcellsf) PetscCall(forest->createcellsf(dm, &forest->cellSF));
1262c7eeac06SToby Isaac   *cellSF = forest->cellSF;
12633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1264c7eeac06SToby Isaac }
1265c7eeac06SToby Isaac 
12669be51f97SToby Isaac /*@C
12679be51f97SToby Isaac   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
1268*a4e35b19SJacob Faibussowitsch   `DMForestGetAdaptivityForest()`) that holds the adaptation flags (refinement, coarsening, or some combination).
12699be51f97SToby Isaac 
127020f4b53cSBarry Smith   Logically Collective
12719be51f97SToby Isaac 
12729be51f97SToby Isaac   Input Parameters:
127360225df5SJacob Faibussowitsch + dm         - the forest
127460225df5SJacob Faibussowitsch - adaptLabel - the label in the pre-adaptation forest
12759be51f97SToby Isaac 
12769be51f97SToby Isaac   Level: intermediate
12779be51f97SToby Isaac 
1278*a4e35b19SJacob Faibussowitsch   Note:
1279*a4e35b19SJacob Faibussowitsch   The interpretation of the label values is up to the subtype of `DMFOREST`, but
1280*a4e35b19SJacob Faibussowitsch   `DM_ADAPT_DETERMINE`, `DM_ADAPT_KEEP`, `DM_ADAPT_REFINE`, and `DM_ADAPT_COARSEN` have been
1281*a4e35b19SJacob Faibussowitsch   reserved as choices that should be accepted by all subtypes.
1282*a4e35b19SJacob Faibussowitsch 
128320f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdaptivityLabel()`
12849be51f97SToby Isaac @*/
1285d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel)
1286d71ae5a4SJacob Faibussowitsch {
1287c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1288c7eeac06SToby Isaac 
1289c7eeac06SToby Isaac   PetscFunctionBegin;
1290c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12911fd50544SStefano Zampini   if (adaptLabel) PetscValidHeaderSpecific(adaptLabel, DMLABEL_CLASSID, 2);
12929566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)adaptLabel));
12939566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&forest->adaptLabel));
12941fd50544SStefano Zampini   forest->adaptLabel = adaptLabel;
12953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1296c7eeac06SToby Isaac }
1297c7eeac06SToby Isaac 
12989be51f97SToby Isaac /*@C
129920f4b53cSBarry Smith   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see `DMForestGetAdaptivityForest()`) that
1300*a4e35b19SJacob Faibussowitsch   holds the adaptation flags (refinement, coarsening, or some combination).
13019be51f97SToby Isaac 
130220f4b53cSBarry Smith   Not Collective
13039be51f97SToby Isaac 
13049be51f97SToby Isaac   Input Parameter:
13059be51f97SToby Isaac . dm - the forest
13069be51f97SToby Isaac 
13079be51f97SToby Isaac   Output Parameter:
13089be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest
13099be51f97SToby Isaac 
13109be51f97SToby Isaac   Level: intermediate
13119be51f97SToby Isaac 
1312*a4e35b19SJacob Faibussowitsch   Note:
1313*a4e35b19SJacob Faibussowitsch   The interpretation of the label values is up to the subtype of `DMFOREST`, but
1314*a4e35b19SJacob Faibussowitsch   `DM_ADAPT_DETERMINE`, `DM_ADAPT_KEEP`, `DM_ADAPT_REFINE`, and `DM_ADAPT_COARSEN` have been
1315*a4e35b19SJacob Faibussowitsch   reserved as choices that should be accepted by all subtypes.
1316*a4e35b19SJacob Faibussowitsch 
131720f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdaptivityLabel()`
13189be51f97SToby Isaac @*/
1319d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel)
1320d71ae5a4SJacob Faibussowitsch {
1321c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1322c7eeac06SToby Isaac 
1323c7eeac06SToby Isaac   PetscFunctionBegin;
1324c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1325ba936b91SToby Isaac   *adaptLabel = forest->adaptLabel;
13263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1327c7eeac06SToby Isaac }
1328c7eeac06SToby Isaac 
13299be51f97SToby Isaac /*@
133020f4b53cSBarry Smith   DMForestSetCellWeights - Set the weights assigned to each of the cells (see `DMForestGetCellChart()`) of the current
133120f4b53cSBarry Smith   process: weights are used to determine parallel partitioning.
13329be51f97SToby Isaac 
133320f4b53cSBarry Smith   Logically Collective
13349be51f97SToby Isaac 
13359be51f97SToby Isaac   Input Parameters:
13369be51f97SToby Isaac + dm       - the forest
133720f4b53cSBarry Smith . weights  - the array of weights (see `DMForestSetWeightCapacity()`) for all cells, or `NULL` to indicate each cell has weight 1.
13389be51f97SToby Isaac - copyMode - how weights should reference weights
13399be51f97SToby Isaac 
13409be51f97SToby Isaac   Level: advanced
13419be51f97SToby Isaac 
134220f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetCellWeights()`, `DMForestSetWeightCapacity()`
13439be51f97SToby Isaac @*/
1344d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1345d71ae5a4SJacob Faibussowitsch {
1346c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1347c7eeac06SToby Isaac   PetscInt   cStart, cEnd;
1348c7eeac06SToby Isaac 
1349c7eeac06SToby Isaac   PetscFunctionBegin;
1350c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
13519566063dSJacob Faibussowitsch   PetscCall(DMForestGetCellChart(dm, &cStart, &cEnd));
135263a3b9bcSJacob Faibussowitsch   PetscCheck(cEnd >= cStart, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "cell chart [%" PetscInt_FMT ",%" PetscInt_FMT ") is not valid", cStart, cEnd);
1353c7eeac06SToby Isaac   if (copyMode == PETSC_COPY_VALUES) {
135448a46eb9SPierre Jolivet     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) PetscCall(PetscMalloc1(cEnd - cStart, &forest->cellWeights));
13559566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(forest->cellWeights, weights, cEnd - cStart));
1356c7eeac06SToby Isaac     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
13573ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1358c7eeac06SToby Isaac   }
135948a46eb9SPierre Jolivet   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) PetscCall(PetscFree(forest->cellWeights));
1360c7eeac06SToby Isaac   forest->cellWeights         = weights;
1361c7eeac06SToby Isaac   forest->cellWeightsCopyMode = copyMode;
13623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1363c7eeac06SToby Isaac }
1364c7eeac06SToby Isaac 
13659be51f97SToby Isaac /*@
136620f4b53cSBarry Smith   DMForestGetCellWeights - Get the weights assigned to each of the cells (see `DMForestGetCellChart()`) of the current
136720f4b53cSBarry Smith   process: weights are used to determine parallel partitioning.
13689be51f97SToby Isaac 
136920f4b53cSBarry Smith   Not Collective
13709be51f97SToby Isaac 
13719be51f97SToby Isaac   Input Parameter:
13729be51f97SToby Isaac . dm - the forest
13739be51f97SToby Isaac 
13749be51f97SToby Isaac   Output Parameter:
137520f4b53cSBarry Smith . weights - the array of weights for all cells, or `NULL` to indicate each cell has weight 1.
13769be51f97SToby Isaac 
13779be51f97SToby Isaac   Level: advanced
13789be51f97SToby Isaac 
137920f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetCellWeights()`, `DMForestSetWeightCapacity()`
13809be51f97SToby Isaac @*/
1381d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1382d71ae5a4SJacob Faibussowitsch {
1383c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1384c7eeac06SToby Isaac 
1385c7eeac06SToby Isaac   PetscFunctionBegin;
1386c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
13874f572ea9SToby Isaac   PetscAssertPointer(weights, 2);
1388c7eeac06SToby Isaac   *weights = forest->cellWeights;
13893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1390c7eeac06SToby Isaac }
1391c7eeac06SToby Isaac 
13929be51f97SToby Isaac /*@
13939be51f97SToby Isaac   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
1394*a4e35b19SJacob Faibussowitsch   a pre-adaptation forest (see `DMForestGetAdaptivityForest()`).
13959be51f97SToby Isaac 
139620f4b53cSBarry Smith   Logically Collective
13979be51f97SToby Isaac 
139860225df5SJacob Faibussowitsch   Input Parameters:
13999be51f97SToby Isaac + dm       - the forest
14009be51f97SToby Isaac - capacity - this process's capacity
14019be51f97SToby Isaac 
14029be51f97SToby Isaac   Level: advanced
14039be51f97SToby Isaac 
1404*a4e35b19SJacob Faibussowitsch   Note:
1405*a4e35b19SJacob Faibussowitsch   After partitioning, the ratio of the weight of each process's cells to the process's capacity
1406*a4e35b19SJacob Faibussowitsch   will be roughly equal for all processes.  A capacity of 0 indicates that the current process
1407*a4e35b19SJacob Faibussowitsch   should not have any cells after repartitioning.
1408*a4e35b19SJacob Faibussowitsch 
1409*a4e35b19SJacob Faibussowitsch .seealso: `DM`, `DMFOREST`, `DMForestGetWeightCapacity()`, `DMForestSetCellWeights()`, `DMForestSetCellWeightFactor()`
14109be51f97SToby Isaac @*/
1411d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1412d71ae5a4SJacob Faibussowitsch {
1413c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1414c7eeac06SToby Isaac 
1415c7eeac06SToby Isaac   PetscFunctionBegin;
1416c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
141728b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the weight capacity after setup");
141863a3b9bcSJacob Faibussowitsch   PetscCheck(capacity >= 0., PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot have negative weight capacity; %g", (double)capacity);
1419c7eeac06SToby Isaac   forest->weightCapacity = capacity;
14203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1421c7eeac06SToby Isaac }
1422c7eeac06SToby Isaac 
14239be51f97SToby Isaac /*@
14249be51f97SToby Isaac   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
1425*a4e35b19SJacob Faibussowitsch   `DMForestGetAdaptivityForest()`).
14269be51f97SToby Isaac 
142720f4b53cSBarry Smith   Not Collective
14289be51f97SToby Isaac 
142960225df5SJacob Faibussowitsch   Input Parameter:
14309be51f97SToby Isaac . dm - the forest
14319be51f97SToby Isaac 
143260225df5SJacob Faibussowitsch   Output Parameter:
14339be51f97SToby Isaac . capacity - this process's capacity
14349be51f97SToby Isaac 
14359be51f97SToby Isaac   Level: advanced
14369be51f97SToby Isaac 
1437*a4e35b19SJacob Faibussowitsch   Note:
1438*a4e35b19SJacob Faibussowitsch   After partitioning, the ratio of the weight of each process's cells to the process's capacity
1439*a4e35b19SJacob Faibussowitsch   will be roughly equal for all processes.  A capacity of 0 indicates that the current process
1440*a4e35b19SJacob Faibussowitsch   should not have any cells after repartitioning.
1441*a4e35b19SJacob Faibussowitsch 
1442*a4e35b19SJacob Faibussowitsch .seealso: `DM`, `DMFOREST`, `DMForestSetWeightCapacity()`, `DMForestSetCellWeights()`, `DMForestSetCellWeightFactor()`
14439be51f97SToby Isaac @*/
1444d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1445d71ae5a4SJacob Faibussowitsch {
1446c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest *)dm->data;
1447c7eeac06SToby Isaac 
1448c7eeac06SToby Isaac   PetscFunctionBegin;
1449c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14504f572ea9SToby Isaac   PetscAssertPointer(capacity, 2);
1451c7eeac06SToby Isaac   *capacity = forest->weightCapacity;
14523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1453c7eeac06SToby Isaac }
1454c7eeac06SToby Isaac 
1455d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(DM dm, PetscOptionItems *PetscOptionsObject)
1456d71ae5a4SJacob Faibussowitsch {
145756ba9f64SToby Isaac   PetscBool                  flg, flg1, flg2, flg3, flg4;
1458dd8e54a2SToby Isaac   DMForestTopology           oldTopo;
1459c7eeac06SToby Isaac   char                       stringBuffer[256];
1460dd8e54a2SToby Isaac   PetscViewer                viewer;
1461dd8e54a2SToby Isaac   PetscViewerFormat          format;
146256ba9f64SToby Isaac   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1463c7eeac06SToby Isaac   PetscReal                  weightsFactor;
1464c7eeac06SToby Isaac   DMForestAdaptivityStrategy adaptStrategy;
1465db4d5e8cSToby Isaac 
1466db4d5e8cSToby Isaac   PetscFunctionBegin;
14679566063dSJacob Faibussowitsch   PetscCall(DMForestGetTopology(dm, &oldTopo));
1468d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "DMForest Options");
14699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-dm_forest_topology", "the topology of the forest's base mesh", "DMForestSetTopology", oldTopo, stringBuffer, sizeof(stringBuffer), &flg1));
14709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsViewer("-dm_forest_base_dm", "load the base DM from a viewer specification", "DMForestSetBaseDM", &viewer, &format, &flg2));
14719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsViewer("-dm_forest_coarse_forest", "load the coarse forest from a viewer specification", "DMForestSetCoarseForest", &viewer, &format, &flg3));
14729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsViewer("-dm_forest_fine_forest", "load the fine forest from a viewer specification", "DMForestSetFineForest", &viewer, &format, &flg4));
14731dca8a05SBarry 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}");
147456ba9f64SToby Isaac   if (flg1) {
14759566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm, (DMForestTopology)stringBuffer));
14769566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm, NULL));
14779566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, NULL));
147856ba9f64SToby Isaac   }
147956ba9f64SToby Isaac   if (flg2) {
1480dd8e54a2SToby Isaac     DM base;
1481dd8e54a2SToby Isaac 
14829566063dSJacob Faibussowitsch     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &base));
14839566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
14849566063dSJacob Faibussowitsch     PetscCall(DMLoad(base, viewer));
14859566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
14869566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm, base));
14879566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&base));
14889566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm, NULL));
14899566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, NULL));
1490dd8e54a2SToby Isaac   }
149156ba9f64SToby Isaac   if (flg3) {
1492dd8e54a2SToby Isaac     DM coarse;
1493dd8e54a2SToby Isaac 
14949566063dSJacob Faibussowitsch     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &coarse));
14959566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
14969566063dSJacob Faibussowitsch     PetscCall(DMLoad(coarse, viewer));
14979566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
14989566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, coarse));
14999566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&coarse));
15009566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm, NULL));
15019566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm, NULL));
1502dd8e54a2SToby Isaac   }
150356ba9f64SToby Isaac   if (flg4) {
1504dd8e54a2SToby Isaac     DM fine;
1505dd8e54a2SToby Isaac 
15069566063dSJacob Faibussowitsch     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &fine));
15079566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
15089566063dSJacob Faibussowitsch     PetscCall(DMLoad(fine, viewer));
15099566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
15109566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm, fine));
15119566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&fine));
15129566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm, NULL));
15139566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm, NULL));
1514dd8e54a2SToby Isaac   }
15159566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdjacencyDimension(dm, &adjDim));
15169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_dimension", "set the dimension of points that define adjacency in the forest", "DMForestSetAdjacencyDimension", adjDim, &adjDim, &flg, 0));
1517dd8e54a2SToby Isaac   if (flg) {
15189566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdjacencyDimension(dm, adjDim));
1519f885a11aSToby Isaac   } else {
15209566063dSJacob Faibussowitsch     PetscCall(DMForestGetAdjacencyCodimension(dm, &adjCodim));
15219566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_codimension", "set the codimension of points that define adjacency in the forest", "DMForestSetAdjacencyCodimension", adjCodim, &adjCodim, &flg, 1));
15221baa6e33SBarry Smith     if (flg) PetscCall(DMForestSetAdjacencyCodimension(dm, adjCodim));
1523dd8e54a2SToby Isaac   }
15249566063dSJacob Faibussowitsch   PetscCall(DMForestGetPartitionOverlap(dm, &overlap));
15259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_partition_overlap", "set the degree of partition overlap", "DMForestSetPartitionOverlap", overlap, &overlap, &flg, 0));
15261baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetPartitionOverlap(dm, overlap));
1527a6121fbdSMatthew G. Knepley #if 0
15289566063dSJacob 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));
1529a6121fbdSMatthew G. Knepley   if (flg) {
15309566063dSJacob Faibussowitsch     PetscCall(DMForestSetMinimumRefinement(dm,minRefinement));
15319566063dSJacob Faibussowitsch     PetscCall(DMForestSetInitialRefinement(dm,minRefinement));
1532a6121fbdSMatthew G. Knepley   }
15339566063dSJacob 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));
1534a6121fbdSMatthew G. Knepley   if (flg) {
15359566063dSJacob Faibussowitsch     PetscCall(DMForestSetMinimumRefinement(dm,0));
15369566063dSJacob Faibussowitsch     PetscCall(DMForestSetInitialRefinement(dm,initRefinement));
1537a6121fbdSMatthew G. Knepley   }
1538a6121fbdSMatthew G. Knepley #endif
15399566063dSJacob Faibussowitsch   PetscCall(DMForestGetMinimumRefinement(dm, &minRefinement));
15409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_minimum_refinement", "set the minimum level of refinement in the forest", "DMForestSetMinimumRefinement", minRefinement, &minRefinement, &flg, 0));
15411baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetMinimumRefinement(dm, minRefinement));
15429566063dSJacob Faibussowitsch   PetscCall(DMForestGetInitialRefinement(dm, &initRefinement));
15439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_initial_refinement", "set the initial level of refinement in the forest", "DMForestSetInitialRefinement", initRefinement, &initRefinement, &flg, 0));
15441baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetInitialRefinement(dm, initRefinement));
15459566063dSJacob Faibussowitsch   PetscCall(DMForestGetMaximumRefinement(dm, &maxRefinement));
15469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_maximum_refinement", "set the maximum level of refinement in the forest", "DMForestSetMaximumRefinement", maxRefinement, &maxRefinement, &flg, 0));
15471baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetMaximumRefinement(dm, maxRefinement));
15489566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivityStrategy(dm, &adaptStrategy));
15499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-dm_forest_adaptivity_strategy", "the forest's adaptivity-flag resolution strategy", "DMForestSetAdaptivityStrategy", adaptStrategy, stringBuffer, sizeof(stringBuffer), &flg));
15501baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetAdaptivityStrategy(dm, (DMForestAdaptivityStrategy)stringBuffer));
15519566063dSJacob Faibussowitsch   PetscCall(DMForestGetGradeFactor(dm, &grade));
15529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_grade_factor", "grade factor between neighboring cells", "DMForestSetGradeFactor", grade, &grade, &flg, 0));
15531baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetGradeFactor(dm, grade));
15549566063dSJacob Faibussowitsch   PetscCall(DMForestGetCellWeightFactor(dm, &weightsFactor));
15559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_forest_cell_weight_factor", "multiplying weight factor for cell refinement", "DMForestSetCellWeightFactor", weightsFactor, &weightsFactor, &flg));
15561baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetCellWeightFactor(dm, weightsFactor));
1557d0609cedSBarry Smith   PetscOptionsHeadEnd();
15583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1559db4d5e8cSToby Isaac }
1560db4d5e8cSToby Isaac 
1561d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1562d71ae5a4SJacob Faibussowitsch {
1563d8984e3bSMatthew G. Knepley   PetscFunctionBegin;
15649566063dSJacob Faibussowitsch   if (subdm) PetscCall(DMClone(dm, subdm));
15659566063dSJacob Faibussowitsch   PetscCall(DMCreateSectionSubDM(dm, numFields, fields, is, subdm));
15663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1567d8984e3bSMatthew G. Knepley }
1568d8984e3bSMatthew G. Knepley 
1569d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
1570d71ae5a4SJacob Faibussowitsch {
15715421bac9SToby Isaac   DMLabel refine;
15725421bac9SToby Isaac   DM      fineDM;
15735421bac9SToby Isaac 
15745421bac9SToby Isaac   PetscFunctionBegin;
15759566063dSJacob Faibussowitsch   PetscCall(DMGetFineDM(dm, &fineDM));
15765421bac9SToby Isaac   if (fineDM) {
15779566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fineDM));
15785421bac9SToby Isaac     *dmRefined = fineDM;
15793ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
15805421bac9SToby Isaac   }
15819566063dSJacob Faibussowitsch   PetscCall(DMForestTemplate(dm, comm, dmRefined));
15829566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "refine", &refine));
15835421bac9SToby Isaac   if (!refine) {
15849566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "refine", &refine));
15859566063dSJacob Faibussowitsch     PetscCall(DMLabelSetDefaultValue(refine, DM_ADAPT_REFINE));
15861baa6e33SBarry Smith   } else PetscCall(PetscObjectReference((PetscObject)refine));
15879566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityLabel(*dmRefined, refine));
15889566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&refine));
15893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15905421bac9SToby Isaac }
15915421bac9SToby Isaac 
1592d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
1593d71ae5a4SJacob Faibussowitsch {
15945421bac9SToby Isaac   DMLabel coarsen;
15955421bac9SToby Isaac   DM      coarseDM;
15965421bac9SToby Isaac 
15975421bac9SToby Isaac   PetscFunctionBegin;
15984098eed7SToby Isaac   {
15994098eed7SToby Isaac     PetscMPIInt mpiComparison;
16004098eed7SToby Isaac     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);
16014098eed7SToby Isaac 
16029566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_compare(comm, dmcomm, &mpiComparison));
16031dca8a05SBarry Smith     PetscCheck(mpiComparison == MPI_IDENT || mpiComparison == MPI_CONGRUENT, dmcomm, PETSC_ERR_SUP, "No support for different communicators yet");
16044098eed7SToby Isaac   }
16059566063dSJacob Faibussowitsch   PetscCall(DMGetCoarseDM(dm, &coarseDM));
16065421bac9SToby Isaac   if (coarseDM) {
16079566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)coarseDM));
16085421bac9SToby Isaac     *dmCoarsened = coarseDM;
16093ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
16105421bac9SToby Isaac   }
16119566063dSJacob Faibussowitsch   PetscCall(DMForestTemplate(dm, comm, dmCoarsened));
16129566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityPurpose(*dmCoarsened, DM_ADAPT_COARSEN));
16139566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "coarsen", &coarsen));
16145421bac9SToby Isaac   if (!coarsen) {
16159566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "coarsen", &coarsen));
16169566063dSJacob Faibussowitsch     PetscCall(DMLabelSetDefaultValue(coarsen, DM_ADAPT_COARSEN));
16171baa6e33SBarry Smith   } else PetscCall(PetscObjectReference((PetscObject)coarsen));
16189566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityLabel(*dmCoarsened, coarsen));
16199566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&coarsen));
16203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16215421bac9SToby Isaac }
16225421bac9SToby Isaac 
1623d71ae5a4SJacob Faibussowitsch PetscErrorCode DMAdaptLabel_Forest(DM dm, PETSC_UNUSED Vec metric, DMLabel label, PETSC_UNUSED DMLabel rgLabel, DM *adaptedDM)
1624d71ae5a4SJacob Faibussowitsch {
162509350103SToby Isaac   PetscBool success;
162609350103SToby Isaac 
162709350103SToby Isaac   PetscFunctionBegin;
16289566063dSJacob Faibussowitsch   PetscCall(DMForestTemplate(dm, PetscObjectComm((PetscObject)dm), adaptedDM));
16299566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityLabel(*adaptedDM, label));
16309566063dSJacob Faibussowitsch   PetscCall(DMSetUp(*adaptedDM));
16319566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivitySuccess(*adaptedDM, &success));
163209350103SToby Isaac   if (!success) {
16339566063dSJacob Faibussowitsch     PetscCall(DMDestroy(adaptedDM));
163409350103SToby Isaac     *adaptedDM = NULL;
163509350103SToby Isaac   }
16363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
163709350103SToby Isaac }
163809350103SToby Isaac 
1639d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMInitialize_Forest(DM dm)
1640d71ae5a4SJacob Faibussowitsch {
1641d222f98bSToby Isaac   PetscFunctionBegin;
16429566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(dm->ops, sizeof(*(dm->ops))));
1643d222f98bSToby Isaac 
1644d222f98bSToby Isaac   dm->ops->clone          = DMClone_Forest;
1645d222f98bSToby Isaac   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1646d222f98bSToby Isaac   dm->ops->destroy        = DMDestroy_Forest;
1647d8984e3bSMatthew G. Knepley   dm->ops->createsubdm    = DMCreateSubDM_Forest;
16485421bac9SToby Isaac   dm->ops->refine         = DMRefine_Forest;
16495421bac9SToby Isaac   dm->ops->coarsen        = DMCoarsen_Forest;
16503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1651d222f98bSToby Isaac }
1652d222f98bSToby Isaac 
16539be51f97SToby Isaac /*MC
16549be51f97SToby Isaac 
165520f4b53cSBarry Smith      DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh.  Forests usually have a base `DM`
165620f4b53cSBarry Smith   (see `DMForestGetBaseDM()`), from which it is refined.  The refinement and partitioning of forests is considered
165720f4b53cSBarry Smith   immutable after `DMSetUp()` is called.  To adapt a mesh, one should call `DMForestTemplate()` to create a new mesh that
165820f4b53cSBarry Smith   will default to being identical to it, specify how that mesh should differ, and then calling `DMSetUp()` on the new
1659bae1f979SBarry Smith   mesh.
1660bae1f979SBarry Smith 
1661bae1f979SBarry Smith   To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the
166220f4b53cSBarry Smith   previous mesh whose values indicate which cells should be refined (`DM_ADAPT_REFINE`) or coarsened (`DM_ADAPT_COARSEN`)
1663bae1f979SBarry Smith   and how (subtypes are free to allow additional values for things like anisotropic refinement).  The label should be
166420f4b53cSBarry Smith   given to the *new* mesh with `DMForestSetAdaptivityLabel()`.
16659be51f97SToby Isaac 
16669be51f97SToby Isaac   Level: advanced
16679be51f97SToby Isaac 
166820f4b53cSBarry Smith .seealso: `DMType`, `DM`, `DMCreate()`, `DMSetType()`, `DMForestGetBaseDM()`, `DMForestSetBaseDM()`, `DMForestTemplate()`, `DMForestSetAdaptivityLabel()`
16699be51f97SToby Isaac M*/
16709be51f97SToby Isaac 
1671d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1672d71ae5a4SJacob Faibussowitsch {
1673db4d5e8cSToby Isaac   DM_Forest *forest;
1674db4d5e8cSToby Isaac 
1675db4d5e8cSToby Isaac   PetscFunctionBegin;
1676db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
16774dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&forest));
1678db4d5e8cSToby Isaac   dm->dim                     = 0;
1679db4d5e8cSToby Isaac   dm->data                    = forest;
1680db4d5e8cSToby Isaac   forest->refct               = 1;
1681db4d5e8cSToby Isaac   forest->data                = NULL;
1682db4d5e8cSToby Isaac   forest->topology            = NULL;
1683cd3c525cSToby Isaac   forest->adapt               = NULL;
1684db4d5e8cSToby Isaac   forest->base                = NULL;
16856a87ffbfSToby Isaac   forest->adaptPurpose        = DM_ADAPT_DETERMINE;
1686db4d5e8cSToby Isaac   forest->adjDim              = PETSC_DEFAULT;
1687db4d5e8cSToby Isaac   forest->overlap             = PETSC_DEFAULT;
1688db4d5e8cSToby Isaac   forest->minRefinement       = PETSC_DEFAULT;
1689db4d5e8cSToby Isaac   forest->maxRefinement       = PETSC_DEFAULT;
169056ba9f64SToby Isaac   forest->initRefinement      = PETSC_DEFAULT;
1691c7eeac06SToby Isaac   forest->cStart              = PETSC_DETERMINE;
1692c7eeac06SToby Isaac   forest->cEnd                = PETSC_DETERMINE;
1693cd3c525cSToby Isaac   forest->cellSF              = NULL;
1694ebdf65a2SToby Isaac   forest->adaptLabel          = NULL;
1695db4d5e8cSToby Isaac   forest->gradeFactor         = 2;
1696db4d5e8cSToby Isaac   forest->cellWeights         = NULL;
1697db4d5e8cSToby Isaac   forest->cellWeightsCopyMode = PETSC_USE_POINTER;
1698db4d5e8cSToby Isaac   forest->weightsFactor       = 1.;
1699db4d5e8cSToby Isaac   forest->weightCapacity      = 1.;
17009566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityStrategy(dm, DMFORESTADAPTALL));
17019566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Forest(dm));
17023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1703db4d5e8cSToby Isaac }
1704