xref: /petsc/src/dm/impls/forest/forest.c (revision 4fb89dddf56594b92bdd2ca7e24874fafe134f45)
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 
1027d4645fSToby Isaac struct _DMForestTypeLink
1127d4645fSToby Isaac {
1227d4645fSToby Isaac   char             *name;
1327d4645fSToby Isaac   DMForestTypeLink next;
1427d4645fSToby Isaac };
1527d4645fSToby Isaac 
1627d4645fSToby Isaac DMForestTypeLink DMForestTypeList;
1727d4645fSToby Isaac 
1827d4645fSToby Isaac static PetscErrorCode DMForestPackageFinalize(void)
1927d4645fSToby Isaac {
2027d4645fSToby Isaac   DMForestTypeLink oldLink, link = DMForestTypeList;
2127d4645fSToby Isaac 
2227d4645fSToby Isaac   PetscFunctionBegin;
2327d4645fSToby Isaac   while (link) {
2427d4645fSToby Isaac     oldLink = link;
259566063dSJacob Faibussowitsch     PetscCall(PetscFree(oldLink->name));
2627d4645fSToby Isaac     link    = oldLink->next;
279566063dSJacob Faibussowitsch     PetscCall(PetscFree(oldLink));
2827d4645fSToby Isaac   }
2927d4645fSToby Isaac   PetscFunctionReturn(0);
3027d4645fSToby Isaac }
3127d4645fSToby Isaac 
3227d4645fSToby Isaac static PetscErrorCode DMForestPackageInitialize(void)
3327d4645fSToby Isaac {
3427d4645fSToby Isaac   PetscFunctionBegin;
3527d4645fSToby Isaac   if (DMForestPackageInitialized) PetscFunctionReturn(0);
3627d4645fSToby Isaac   DMForestPackageInitialized = PETSC_TRUE;
37f885a11aSToby Isaac 
389566063dSJacob Faibussowitsch   PetscCall(DMForestRegisterType(DMFOREST));
399566063dSJacob Faibussowitsch   PetscCall(PetscRegisterFinalize(DMForestPackageFinalize));
4027d4645fSToby Isaac   PetscFunctionReturn(0);
4127d4645fSToby Isaac }
4227d4645fSToby Isaac 
439be51f97SToby Isaac /*@C
449be51f97SToby Isaac   DMForestRegisterType - Registers a DMType as a subtype of DMFOREST (so that DMIsForest() will be correct)
459be51f97SToby Isaac 
469be51f97SToby Isaac   Not Collective
479be51f97SToby Isaac 
489be51f97SToby Isaac   Input parameter:
499be51f97SToby Isaac . name - the name of the type
509be51f97SToby Isaac 
519be51f97SToby Isaac   Level: advanced
529be51f97SToby Isaac 
53db781477SPatrick Sanan .seealso: `DMFOREST`, `DMIsForest()`
549be51f97SToby Isaac @*/
5527d4645fSToby Isaac PetscErrorCode DMForestRegisterType(DMType name)
5627d4645fSToby Isaac {
5727d4645fSToby Isaac   DMForestTypeLink link;
5827d4645fSToby Isaac 
5927d4645fSToby Isaac   PetscFunctionBegin;
609566063dSJacob Faibussowitsch   PetscCall(DMForestPackageInitialize());
619566063dSJacob Faibussowitsch   PetscCall(PetscNew(&link));
629566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name,&link->name));
6327d4645fSToby Isaac   link->next       = DMForestTypeList;
6427d4645fSToby Isaac   DMForestTypeList = link;
6527d4645fSToby Isaac   PetscFunctionReturn(0);
6627d4645fSToby Isaac }
6727d4645fSToby Isaac 
689be51f97SToby Isaac /*@
699be51f97SToby Isaac   DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes
709be51f97SToby Isaac 
719be51f97SToby Isaac   Not Collective
729be51f97SToby Isaac 
739be51f97SToby Isaac   Input parameter:
749be51f97SToby Isaac . dm - the DM object
759be51f97SToby Isaac 
769be51f97SToby Isaac   Output parameter:
779be51f97SToby Isaac . isForest - whether dm is a subtype of DMFOREST
789be51f97SToby Isaac 
799be51f97SToby Isaac   Level: intermediate
809be51f97SToby Isaac 
81db781477SPatrick Sanan .seealso: `DMFOREST`, `DMForestRegisterType()`
829be51f97SToby Isaac @*/
8327d4645fSToby Isaac PetscErrorCode DMIsForest(DM dm, PetscBool *isForest)
8427d4645fSToby Isaac {
8527d4645fSToby Isaac   DMForestTypeLink link = DMForestTypeList;
8627d4645fSToby Isaac 
8727d4645fSToby Isaac   PetscFunctionBegin;
8827d4645fSToby Isaac   while (link) {
8927d4645fSToby Isaac     PetscBool sameType;
909566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)dm,link->name,&sameType));
9127d4645fSToby Isaac     if (sameType) {
9227d4645fSToby Isaac       *isForest = PETSC_TRUE;
9327d4645fSToby Isaac       PetscFunctionReturn(0);
9427d4645fSToby Isaac     }
9527d4645fSToby Isaac     link = link->next;
9627d4645fSToby Isaac   }
9727d4645fSToby Isaac   *isForest = PETSC_FALSE;
9827d4645fSToby Isaac   PetscFunctionReturn(0);
9927d4645fSToby Isaac }
10027d4645fSToby Isaac 
1019be51f97SToby Isaac /*@
1029be51f97SToby Isaac   DMForestTemplate - Create a new DM that will be adapted from a source DM.  The new DM reproduces the configuration
1039be51f97SToby Isaac   of the source, but is not yet setup, so that the user can then define only the ways that the new DM should differ
1049be51f97SToby Isaac   (by, e.g., refinement or repartitioning).  The source DM is also set as the adaptivity source DM of the new DM (see
1059be51f97SToby Isaac   DMForestSetAdaptivityForest()).
1069be51f97SToby Isaac 
1079be51f97SToby Isaac   Collective on dm
1089be51f97SToby Isaac 
1099be51f97SToby Isaac   Input Parameters:
1109be51f97SToby Isaac + dm - the source DM object
1119be51f97SToby Isaac - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen())
1129be51f97SToby Isaac 
1139be51f97SToby Isaac   Output Parameter:
1149be51f97SToby Isaac . tdm - the new DM object
1159be51f97SToby Isaac 
1169be51f97SToby Isaac   Level: intermediate
1179be51f97SToby Isaac 
118db781477SPatrick Sanan .seealso: `DMForestSetAdaptivityForest()`
1199be51f97SToby Isaac @*/
1209be51f97SToby Isaac PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm)
121a0452a8eSToby Isaac {
122a0452a8eSToby Isaac   DM_Forest                  *forest = (DM_Forest*) dm->data;
12320e8089bSToby Isaac   DMType                     type;
124a0452a8eSToby Isaac   DM                         base;
125a0452a8eSToby Isaac   DMForestTopology           topology;
12605e99e11SStefano Zampini   MatType                    mtype;
127a0452a8eSToby Isaac   PetscInt                   dim, overlap, ref, factor;
128a0452a8eSToby Isaac   DMForestAdaptivityStrategy strat;
129795844e7SToby Isaac   void                       *ctx;
13049fc9a2fSToby Isaac   PetscErrorCode             (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void*);
1313e58adeeSToby Isaac   void                       *mapCtx;
132a0452a8eSToby Isaac 
133a0452a8eSToby Isaac   PetscFunctionBegin;
134a0452a8eSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1359566063dSJacob Faibussowitsch   PetscCall(DMCreate(PetscObjectComm((PetscObject)dm),tdm));
1369566063dSJacob Faibussowitsch   PetscCall(DMGetType(dm,&type));
1379566063dSJacob Faibussowitsch   PetscCall(DMSetType(*tdm,type));
1389566063dSJacob Faibussowitsch   PetscCall(DMForestGetBaseDM(dm,&base));
1399566063dSJacob Faibussowitsch   PetscCall(DMForestSetBaseDM(*tdm,base));
1409566063dSJacob Faibussowitsch   PetscCall(DMForestGetTopology(dm,&topology));
1419566063dSJacob Faibussowitsch   PetscCall(DMForestSetTopology(*tdm,topology));
1429566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdjacencyDimension(dm,&dim));
1439566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdjacencyDimension(*tdm,dim));
1449566063dSJacob Faibussowitsch   PetscCall(DMForestGetPartitionOverlap(dm,&overlap));
1459566063dSJacob Faibussowitsch   PetscCall(DMForestSetPartitionOverlap(*tdm,overlap));
1469566063dSJacob Faibussowitsch   PetscCall(DMForestGetMinimumRefinement(dm,&ref));
1479566063dSJacob Faibussowitsch   PetscCall(DMForestSetMinimumRefinement(*tdm,ref));
1489566063dSJacob Faibussowitsch   PetscCall(DMForestGetMaximumRefinement(dm,&ref));
1499566063dSJacob Faibussowitsch   PetscCall(DMForestSetMaximumRefinement(*tdm,ref));
1509566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivityStrategy(dm,&strat));
1519566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityStrategy(*tdm,strat));
1529566063dSJacob Faibussowitsch   PetscCall(DMForestGetGradeFactor(dm,&factor));
1539566063dSJacob Faibussowitsch   PetscCall(DMForestSetGradeFactor(*tdm,factor));
1549566063dSJacob Faibussowitsch   PetscCall(DMForestGetBaseCoordinateMapping(dm,&map,&mapCtx));
1559566063dSJacob Faibussowitsch   PetscCall(DMForestSetBaseCoordinateMapping(*tdm,map,mapCtx));
1561baa6e33SBarry Smith   if (forest->ftemplate) PetscCall((*forest->ftemplate)(dm, *tdm));
1579566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityForest(*tdm,dm));
1589566063dSJacob Faibussowitsch   PetscCall(DMCopyDisc(dm,*tdm));
1599566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(dm,&ctx));
1609566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(*tdm,&ctx));
16190b157c4SStefano Zampini   {
162*4fb89dddSMatthew G. Knepley     const PetscReal *maxCell, *L, *Lstart;
163795844e7SToby Isaac 
164*4fb89dddSMatthew G. Knepley     PetscCall(DMGetPeriodicity(dm,  &maxCell, &Lstart, &L));
165*4fb89dddSMatthew G. Knepley     PetscCall(DMSetPeriodicity(*tdm, maxCell,  Lstart,  L));
166795844e7SToby Isaac   }
1679566063dSJacob Faibussowitsch   PetscCall(DMGetMatType(dm,&mtype));
1689566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(*tdm,mtype));
169a0452a8eSToby Isaac   PetscFunctionReturn(0);
170a0452a8eSToby Isaac }
171a0452a8eSToby Isaac 
17201d9d024SToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm);
17301d9d024SToby Isaac 
174db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm)
175db4d5e8cSToby Isaac {
176db4d5e8cSToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
177db4d5e8cSToby Isaac   const char     *type;
178db4d5e8cSToby Isaac 
179db4d5e8cSToby Isaac   PetscFunctionBegin;
180db4d5e8cSToby Isaac   forest->refct++;
181db4d5e8cSToby Isaac   (*newdm)->data = forest;
1829566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetType((PetscObject) dm, &type));
1839566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject) *newdm, type));
1849566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Forest(*newdm));
185db4d5e8cSToby Isaac   PetscFunctionReturn(0);
186db4d5e8cSToby Isaac }
187db4d5e8cSToby Isaac 
188d222f98bSToby Isaac static PetscErrorCode DMDestroy_Forest(DM dm)
189db4d5e8cSToby Isaac {
190db4d5e8cSToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
191db4d5e8cSToby Isaac 
192db4d5e8cSToby Isaac   PetscFunctionBegin;
193db4d5e8cSToby Isaac   if (--forest->refct > 0) PetscFunctionReturn(0);
1949566063dSJacob Faibussowitsch   if (forest->destroy) PetscCall((*forest->destroy)(dm));
1959566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&forest->cellSF));
1969566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&forest->preCoarseToFine));
1979566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&forest->coarseToPreFine));
1989566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&forest->adaptLabel));
1999566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest->adaptStrategy));
2009566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&forest->base));
2019566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&forest->adapt));
2029566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest->topology));
2039566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest));
204db4d5e8cSToby Isaac   PetscFunctionReturn(0);
205db4d5e8cSToby Isaac }
206db4d5e8cSToby Isaac 
2079be51f97SToby Isaac /*@C
2089be51f97SToby Isaac   DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase.  The topology is a string (e.g.
20968d54884SBarry Smith   "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest during
2109be51f97SToby Isaac   DMSetUp().
2119be51f97SToby Isaac 
2129be51f97SToby Isaac   Logically collective on dm
2139be51f97SToby Isaac 
2149be51f97SToby Isaac   Input parameters:
2159be51f97SToby Isaac + dm - the forest
2169be51f97SToby Isaac - topology - the topology of the forest
2179be51f97SToby Isaac 
2189be51f97SToby Isaac   Level: intermediate
2199be51f97SToby Isaac 
220db781477SPatrick Sanan .seealso: `DMForestGetTopology()`, `DMForestSetBaseDM()`
2219be51f97SToby Isaac @*/
222dd8e54a2SToby Isaac PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology)
223db4d5e8cSToby Isaac {
224db4d5e8cSToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
225db4d5e8cSToby Isaac 
226db4d5e8cSToby Isaac   PetscFunctionBegin;
227db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
22828b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup");
2299566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest->topology));
2309566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy((const char*)topology,(char**) &forest->topology));
231db4d5e8cSToby Isaac   PetscFunctionReturn(0);
232db4d5e8cSToby Isaac }
233db4d5e8cSToby Isaac 
2349be51f97SToby Isaac /*@C
2359be51f97SToby Isaac   DMForestGetTopology - Get a string describing the topology of a DMForest.
2369be51f97SToby Isaac 
2379be51f97SToby Isaac   Not collective
2389be51f97SToby Isaac 
2399be51f97SToby Isaac   Input parameter:
2409be51f97SToby Isaac . dm - the forest
2419be51f97SToby Isaac 
2429be51f97SToby Isaac   Output parameter:
2439be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell')
2449be51f97SToby Isaac 
2459be51f97SToby Isaac   Level: intermediate
2469be51f97SToby Isaac 
247db781477SPatrick Sanan .seealso: `DMForestSetTopology()`
2489be51f97SToby Isaac @*/
249dd8e54a2SToby Isaac PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology)
250dd8e54a2SToby Isaac {
251dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
252dd8e54a2SToby Isaac 
253dd8e54a2SToby Isaac   PetscFunctionBegin;
254dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
255dd8e54a2SToby Isaac   PetscValidPointer(topology,2);
256dd8e54a2SToby Isaac   *topology = forest->topology;
257dd8e54a2SToby Isaac   PetscFunctionReturn(0);
258dd8e54a2SToby Isaac }
259dd8e54a2SToby Isaac 
2609be51f97SToby Isaac /*@
2619be51f97SToby Isaac   DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest.  The
2629be51f97SToby Isaac   forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its
263765b024eSBarry Smith   base.  In general, two forest must share a base to be comparable, to do things like construct interpolators.
2649be51f97SToby Isaac 
2659be51f97SToby Isaac   Logically collective on dm
2669be51f97SToby Isaac 
2679be51f97SToby Isaac   Input Parameters:
2689be51f97SToby Isaac + dm - the forest
2699be51f97SToby Isaac - base - the base DM of the forest
2709be51f97SToby Isaac 
27195452b02SPatrick Sanan   Notes:
27295452b02SPatrick Sanan     Currently the base DM must be a DMPLEX
273765b024eSBarry Smith 
2749be51f97SToby Isaac   Level: intermediate
2759be51f97SToby Isaac 
276db781477SPatrick Sanan .seealso: `DMForestGetBaseDM()`
2779be51f97SToby Isaac @*/
278dd8e54a2SToby Isaac PetscErrorCode DMForestSetBaseDM(DM dm, DM base)
279dd8e54a2SToby Isaac {
280dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
281dd8e54a2SToby Isaac   PetscInt       dim, dimEmbed;
282dd8e54a2SToby Isaac 
283dd8e54a2SToby Isaac   PetscFunctionBegin;
284dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
28528b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup");
2869566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)base));
2879566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&forest->base));
288dd8e54a2SToby Isaac   forest->base = base;
289a0452a8eSToby Isaac   if (base) {
290*4fb89dddSMatthew G. Knepley     const PetscReal *maxCell, *Lstart, *L;
29128dfcf7cSStefano Zampini 
292a0452a8eSToby Isaac     PetscValidHeaderSpecific(base, DM_CLASSID, 2);
2939566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(base,&dim));
2949566063dSJacob Faibussowitsch     PetscCall(DMSetDimension(dm,dim));
2959566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDim(base,&dimEmbed));
2969566063dSJacob Faibussowitsch     PetscCall(DMSetCoordinateDim(dm,dimEmbed));
297*4fb89dddSMatthew G. Knepley     PetscCall(DMGetPeriodicity(base,&maxCell,&Lstart,&L));
298*4fb89dddSMatthew G. Knepley     PetscCall(DMSetPeriodicity(dm,maxCell,Lstart,L));
299*4fb89dddSMatthew G. Knepley   } else PetscCall(DMSetPeriodicity(dm,NULL,NULL,NULL));
300dd8e54a2SToby Isaac   PetscFunctionReturn(0);
301dd8e54a2SToby Isaac }
302dd8e54a2SToby Isaac 
3039be51f97SToby Isaac /*@
3049be51f97SToby Isaac   DMForestGetBaseDM - Get the base DM of a DMForest forest.  The forest will be hierarchically refined from the base,
30568d54884SBarry Smith   and all refinements/coarsenings of the forest will share its base.  In general, two forest must share a base to be
3069be51f97SToby Isaac   comparable, to do things like construct interpolators.
3079be51f97SToby Isaac 
3089be51f97SToby Isaac   Not collective
3099be51f97SToby Isaac 
3109be51f97SToby Isaac   Input Parameter:
3119be51f97SToby Isaac . dm - the forest
3129be51f97SToby Isaac 
3139be51f97SToby Isaac   Output Parameter:
3149be51f97SToby Isaac . base - the base DM of the forest
3159be51f97SToby Isaac 
316367003a6SStefano Zampini   Notes:
317367003a6SStefano Zampini     After DMSetUp(), the base DM will be redundantly distributed across MPI processes
318367003a6SStefano Zampini 
3199be51f97SToby Isaac   Level: intermediate
3209be51f97SToby Isaac 
321db781477SPatrick Sanan .seealso: `DMForestSetBaseDM()`
3229be51f97SToby Isaac @*/
323dd8e54a2SToby Isaac PetscErrorCode DMForestGetBaseDM(DM dm, DM *base)
324dd8e54a2SToby Isaac {
325dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
326dd8e54a2SToby Isaac 
327dd8e54a2SToby Isaac   PetscFunctionBegin;
328dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
329dd8e54a2SToby Isaac   PetscValidPointer(base, 2);
330dd8e54a2SToby Isaac   *base = forest->base;
331dd8e54a2SToby Isaac   PetscFunctionReturn(0);
332dd8e54a2SToby Isaac }
333dd8e54a2SToby Isaac 
33499478f86SToby Isaac PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
335cf38a08cSToby Isaac {
336cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
337cf38a08cSToby Isaac 
338cf38a08cSToby Isaac   PetscFunctionBegin;
339cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
340cf38a08cSToby Isaac   forest->mapcoordinates    = func;
341cf38a08cSToby Isaac   forest->mapcoordinatesctx = ctx;
342cf38a08cSToby Isaac   PetscFunctionReturn(0);
343cf38a08cSToby Isaac }
344cf38a08cSToby Isaac 
34599478f86SToby Isaac PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func) (DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
346cf38a08cSToby Isaac {
347cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
348cf38a08cSToby Isaac 
349cf38a08cSToby Isaac   PetscFunctionBegin;
350cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
351cf38a08cSToby Isaac   if (func) *func = forest->mapcoordinates;
352cf38a08cSToby Isaac   if (ctx) *((void**) ctx) = forest->mapcoordinatesctx;
353cf38a08cSToby Isaac   PetscFunctionReturn(0);
354cf38a08cSToby Isaac }
355cf38a08cSToby Isaac 
3569be51f97SToby Isaac /*@
3579be51f97SToby Isaac   DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be
3589be51f97SToby Isaac   adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp().  Usually not needed
3599be51f97SToby Isaac   by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this
3609be51f97SToby Isaac   routine.
3619be51f97SToby Isaac 
362dffe73a3SToby Isaac   Note that this can be called after setup with adapt = NULL, which will clear all internal data related to the
363dffe73a3SToby Isaac   adaptivity forest from dm.  This way, repeatedly adapting does not leave stale DM objects in memory.
364dffe73a3SToby Isaac 
3659be51f97SToby Isaac   Logically collective on dm
3669be51f97SToby Isaac 
367d8d19677SJose E. Roman   Input Parameters:
3689be51f97SToby Isaac + dm - the new forest, which will be constructed from adapt
3699be51f97SToby Isaac - adapt - the old forest
3709be51f97SToby Isaac 
3719be51f97SToby Isaac   Level: intermediate
3729be51f97SToby Isaac 
373db781477SPatrick Sanan .seealso: `DMForestGetAdaptivityForest()`, `DMForestSetAdaptivityPurpose()`
3749be51f97SToby Isaac @*/
375ba936b91SToby Isaac PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt)
376dd8e54a2SToby Isaac {
377dffe73a3SToby Isaac   DM_Forest      *forest, *adaptForest, *oldAdaptForest;
378dffe73a3SToby Isaac   DM             oldAdapt;
379456cc5b7SMatthew G. Knepley   PetscBool      isForest;
380dd8e54a2SToby Isaac 
381dd8e54a2SToby Isaac   PetscFunctionBegin;
382dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3831fd50544SStefano Zampini   if (adapt) PetscValidHeaderSpecific(adapt, DM_CLASSID, 2);
3849566063dSJacob Faibussowitsch   PetscCall(DMIsForest(dm, &isForest));
385456cc5b7SMatthew G. Knepley   if (!isForest) PetscFunctionReturn(0);
3861dca8a05SBarry Smith   PetscCheck(adapt == NULL || !dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
387ba936b91SToby Isaac   forest         = (DM_Forest*) dm->data;
3889566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivityForest(dm,&oldAdapt));
389193eb951SToby Isaac   adaptForest    = (DM_Forest*) (adapt ? adapt->data : NULL);
390193eb951SToby Isaac   oldAdaptForest = (DM_Forest*) (oldAdapt ? oldAdapt->data : NULL);
391dffe73a3SToby Isaac   if (adaptForest != oldAdaptForest) {
3929566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&forest->preCoarseToFine));
3939566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&forest->coarseToPreFine));
3949566063dSJacob Faibussowitsch     if (forest->clearadaptivityforest) PetscCall((*forest->clearadaptivityforest)(dm));
395dffe73a3SToby Isaac   }
39626d9498aSToby Isaac   switch (forest->adaptPurpose) {
397cd3c525cSToby Isaac   case DM_ADAPT_DETERMINE:
3989566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)adapt));
3999566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&(forest->adapt)));
400ba936b91SToby Isaac     forest->adapt = adapt;
40126d9498aSToby Isaac     break;
402a1b0c543SToby Isaac   case DM_ADAPT_REFINE:
4039566063dSJacob Faibussowitsch     PetscCall(DMSetCoarseDM(dm,adapt));
40426d9498aSToby Isaac     break;
405a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
406bf2d5fbbSStefano Zampini   case DM_ADAPT_COARSEN_LAST:
4079566063dSJacob Faibussowitsch     PetscCall(DMSetFineDM(dm,adapt));
40826d9498aSToby Isaac     break;
40926d9498aSToby Isaac   default:
41026d9498aSToby Isaac     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
41126d9498aSToby Isaac   }
412dd8e54a2SToby Isaac   PetscFunctionReturn(0);
413dd8e54a2SToby Isaac }
414dd8e54a2SToby Isaac 
4159be51f97SToby Isaac /*@
4169be51f97SToby Isaac   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.
4179be51f97SToby Isaac 
4189be51f97SToby Isaac   Not collective
4199be51f97SToby Isaac 
4209be51f97SToby Isaac   Input Parameter:
4219be51f97SToby Isaac . dm - the forest
4229be51f97SToby Isaac 
4239be51f97SToby Isaac   Output Parameter:
4249be51f97SToby Isaac . adapt - the forest from which dm is/was adapted
4259be51f97SToby Isaac 
4269be51f97SToby Isaac   Level: intermediate
4279be51f97SToby Isaac 
428db781477SPatrick Sanan .seealso: `DMForestSetAdaptivityForest()`, `DMForestSetAdaptivityPurpose()`
4299be51f97SToby Isaac @*/
430ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
431dd8e54a2SToby Isaac {
432ba936b91SToby Isaac   DM_Forest      *forest;
433dd8e54a2SToby Isaac 
434dd8e54a2SToby Isaac   PetscFunctionBegin;
435dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
436ba936b91SToby Isaac   forest = (DM_Forest*) dm->data;
43726d9498aSToby Isaac   switch (forest->adaptPurpose) {
438cd3c525cSToby Isaac   case DM_ADAPT_DETERMINE:
439ba936b91SToby Isaac     *adapt = forest->adapt;
44026d9498aSToby Isaac     break;
441a1b0c543SToby Isaac   case DM_ADAPT_REFINE:
4429566063dSJacob Faibussowitsch     PetscCall(DMGetCoarseDM(dm,adapt));
44326d9498aSToby Isaac     break;
444a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
445bf2d5fbbSStefano Zampini   case DM_ADAPT_COARSEN_LAST:
4469566063dSJacob Faibussowitsch     PetscCall(DMGetFineDM(dm,adapt));
44726d9498aSToby Isaac     break;
44826d9498aSToby Isaac   default:
44926d9498aSToby Isaac     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
45026d9498aSToby Isaac   }
45126d9498aSToby Isaac   PetscFunctionReturn(0);
45226d9498aSToby Isaac }
45326d9498aSToby Isaac 
4549be51f97SToby Isaac /*@
4559be51f97SToby Isaac   DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its
456a1b0c543SToby Isaac   source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening
457cd3c525cSToby Isaac   (DM_ADAPT_COARSEN), or undefined (DM_ADAPT_DETERMINE).  This only matters for the purposes of reference counting:
4589be51f97SToby Isaac   during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse
4599be51f97SToby Isaac   relationship (see DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does
4609be51f97SToby Isaac   not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can
4619be51f97SToby Isaac   cause a memory leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
4629be51f97SToby Isaac 
4639be51f97SToby Isaac   Logically collective on dm
4649be51f97SToby Isaac 
4659be51f97SToby Isaac   Input Parameters:
4669be51f97SToby Isaac + dm - the forest
467bf2d5fbbSStefano Zampini - purpose - the adaptivity purpose
4689be51f97SToby Isaac 
4699be51f97SToby Isaac   Level: advanced
4709be51f97SToby Isaac 
471db781477SPatrick Sanan .seealso: `DMForestTemplate()`, `DMForestSetAdaptivityForest()`, `DMForestGetAdaptivityForest()`, `DMAdaptFlag`
4729be51f97SToby Isaac @*/
473a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose)
47426d9498aSToby Isaac {
47526d9498aSToby Isaac   DM_Forest      *forest;
47626d9498aSToby Isaac 
47726d9498aSToby Isaac   PetscFunctionBegin;
47826d9498aSToby Isaac   forest = (DM_Forest*) dm->data;
47928b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
48026d9498aSToby Isaac   if (purpose != forest->adaptPurpose) {
48126d9498aSToby Isaac     DM adapt;
48226d9498aSToby Isaac 
4839566063dSJacob Faibussowitsch     PetscCall(DMForestGetAdaptivityForest(dm,&adapt));
4849566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)adapt));
4859566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm,NULL));
486f885a11aSToby Isaac 
48726d9498aSToby Isaac     forest->adaptPurpose = purpose;
488f885a11aSToby Isaac 
4899566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm,adapt));
4909566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&adapt));
49126d9498aSToby Isaac   }
492dd8e54a2SToby Isaac   PetscFunctionReturn(0);
493dd8e54a2SToby Isaac }
494dd8e54a2SToby Isaac 
49556c0450aSToby Isaac /*@
49656c0450aSToby Isaac   DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with
497bf2d5fbbSStefano Zampini   DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening (DM_ADAPT_COARSEN),
498bf2d5fbbSStefano Zampini   coarsening only the last level (DM_ADAPT_COARSEN_LAST) or undefined (DM_ADAPT_DETERMINE).
499bf2d5fbbSStefano Zampini   This only matters for the purposes of reference counting: during DMDestroy(), cyclic
50056c0450aSToby Isaac   references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see
50156c0450aSToby Isaac   DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does not maintain a
50256c0450aSToby Isaac   reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory
50356c0450aSToby Isaac   leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
50456c0450aSToby Isaac 
50556c0450aSToby Isaac   Not collective
50656c0450aSToby Isaac 
50756c0450aSToby Isaac   Input Parameter:
50856c0450aSToby Isaac . dm - the forest
50956c0450aSToby Isaac 
51056c0450aSToby Isaac   Output Parameter:
511bf2d5fbbSStefano Zampini . purpose - the adaptivity purpose
51256c0450aSToby Isaac 
51356c0450aSToby Isaac   Level: advanced
51456c0450aSToby Isaac 
515db781477SPatrick Sanan .seealso: `DMForestTemplate()`, `DMForestSetAdaptivityForest()`, `DMForestGetAdaptivityForest()`, `DMAdaptFlag`
51656c0450aSToby Isaac @*/
517a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose)
51856c0450aSToby Isaac {
51956c0450aSToby Isaac   DM_Forest *forest;
52056c0450aSToby Isaac 
52156c0450aSToby Isaac   PetscFunctionBegin;
52256c0450aSToby Isaac   forest   = (DM_Forest*) dm->data;
52356c0450aSToby Isaac   *purpose = forest->adaptPurpose;
52456c0450aSToby Isaac   PetscFunctionReturn(0);
52556c0450aSToby Isaac }
52656c0450aSToby Isaac 
5279be51f97SToby Isaac /*@
5289be51f97SToby Isaac   DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
5299be51f97SToby Isaac   cell adjacency (for the purposes of partitioning and overlap).
5309be51f97SToby Isaac 
5319be51f97SToby Isaac   Logically collective on dm
5329be51f97SToby Isaac 
5339be51f97SToby Isaac   Input Parameters:
5349be51f97SToby Isaac + dm - the forest
5359be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency)
5369be51f97SToby Isaac 
5379be51f97SToby Isaac   Level: intermediate
5389be51f97SToby Isaac 
539db781477SPatrick Sanan .seealso: `DMForestGetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`, `DMForestSetPartitionOverlap()`
5409be51f97SToby Isaac @*/
541dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
542dd8e54a2SToby Isaac {
543dd8e54a2SToby Isaac   PetscInt       dim;
544dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
545dd8e54a2SToby Isaac 
546dd8e54a2SToby Isaac   PetscFunctionBegin;
547dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
54828b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup");
54963a3b9bcSJacob Faibussowitsch   PetscCheck(adjDim >= 0,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %" PetscInt_FMT, adjDim);
5509566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
55163a3b9bcSJacob Faibussowitsch   PetscCheck(adjDim <= dim,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %" PetscInt_FMT ": %" PetscInt_FMT, dim, adjDim);
552dd8e54a2SToby Isaac   forest->adjDim = adjDim;
553dd8e54a2SToby Isaac   PetscFunctionReturn(0);
554dd8e54a2SToby Isaac }
555dd8e54a2SToby Isaac 
5569be51f97SToby Isaac /*@
5579be51f97SToby Isaac   DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that,
5589be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
5599be51f97SToby Isaac 
5609be51f97SToby Isaac   Logically collective on dm
5619be51f97SToby Isaac 
5629be51f97SToby Isaac   Input Parameters:
5639be51f97SToby Isaac + dm - the forest
5649be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
5659be51f97SToby Isaac 
5669be51f97SToby Isaac   Level: intermediate
5679be51f97SToby Isaac 
568db781477SPatrick Sanan .seealso: `DMForestGetAdjacencyCodimension()`, `DMForestSetAdjacencyDimension()`
5699be51f97SToby Isaac @*/
570dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
571dd8e54a2SToby Isaac {
572dd8e54a2SToby Isaac   PetscInt       dim;
573dd8e54a2SToby Isaac 
574dd8e54a2SToby Isaac   PetscFunctionBegin;
575dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5769566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
5779566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdjacencyDimension(dm,dim-adjCodim));
578dd8e54a2SToby Isaac   PetscFunctionReturn(0);
579dd8e54a2SToby Isaac }
580dd8e54a2SToby Isaac 
5819be51f97SToby Isaac /*@
5829be51f97SToby Isaac   DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
5839be51f97SToby Isaac   purposes of partitioning and overlap).
5849be51f97SToby Isaac 
5859be51f97SToby Isaac   Not collective
5869be51f97SToby Isaac 
5879be51f97SToby Isaac   Input Parameter:
5889be51f97SToby Isaac . dm - the forest
5899be51f97SToby Isaac 
5909be51f97SToby Isaac   Output Parameter:
5919be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency)
5929be51f97SToby Isaac 
5939be51f97SToby Isaac   Level: intermediate
5949be51f97SToby Isaac 
595db781477SPatrick Sanan .seealso: `DMForestSetAdjacencyDimension()`, `DMForestGetAdjacencyCodimension()`, `DMForestSetPartitionOverlap()`
5969be51f97SToby Isaac @*/
597dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
598dd8e54a2SToby Isaac {
599dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
600dd8e54a2SToby Isaac 
601dd8e54a2SToby Isaac   PetscFunctionBegin;
602dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
603dd8e54a2SToby Isaac   PetscValidIntPointer(adjDim,2);
604dd8e54a2SToby Isaac   *adjDim = forest->adjDim;
605dd8e54a2SToby Isaac   PetscFunctionReturn(0);
606dd8e54a2SToby Isaac }
607dd8e54a2SToby Isaac 
6089be51f97SToby Isaac /*@
6099be51f97SToby Isaac   DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that,
6109be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
6119be51f97SToby Isaac 
6129be51f97SToby Isaac   Not collective
6139be51f97SToby Isaac 
6149be51f97SToby Isaac   Input Parameter:
6159be51f97SToby Isaac . dm - the forest
6169be51f97SToby Isaac 
6179be51f97SToby Isaac   Output Parameter:
6189be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
6199be51f97SToby Isaac 
6209be51f97SToby Isaac   Level: intermediate
6219be51f97SToby Isaac 
622db781477SPatrick Sanan .seealso: `DMForestSetAdjacencyCodimension()`, `DMForestGetAdjacencyDimension()`
6239be51f97SToby Isaac @*/
624dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
625dd8e54a2SToby Isaac {
626dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
627dd8e54a2SToby Isaac   PetscInt       dim;
628dd8e54a2SToby Isaac 
629dd8e54a2SToby Isaac   PetscFunctionBegin;
630dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
631dd8e54a2SToby Isaac   PetscValidIntPointer(adjCodim,2);
6329566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm,&dim));
633dd8e54a2SToby Isaac   *adjCodim = dim - forest->adjDim;
634dd8e54a2SToby Isaac   PetscFunctionReturn(0);
635dd8e54a2SToby Isaac }
636dd8e54a2SToby Isaac 
6379be51f97SToby Isaac /*@
6389be51f97SToby Isaac   DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel
6399be51f97SToby Isaac   partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding
6409be51f97SToby Isaac   adjacent cells
6419be51f97SToby Isaac 
6429be51f97SToby Isaac   Logically collective on dm
6439be51f97SToby Isaac 
6449be51f97SToby Isaac   Input Parameters:
6459be51f97SToby Isaac + dm - the forest
6469be51f97SToby Isaac - overlap - default 0
6479be51f97SToby Isaac 
6489be51f97SToby Isaac   Level: intermediate
6499be51f97SToby Isaac 
650db781477SPatrick Sanan .seealso: `DMForestGetPartitionOverlap()`, `DMForestSetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`
6519be51f97SToby Isaac @*/
652dd8e54a2SToby Isaac PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap)
653dd8e54a2SToby Isaac {
654dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
655dd8e54a2SToby Isaac 
656dd8e54a2SToby Isaac   PetscFunctionBegin;
657dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
65828b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup");
65963a3b9bcSJacob Faibussowitsch   PetscCheck(overlap >= 0,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %" PetscInt_FMT, overlap);
660dd8e54a2SToby Isaac   forest->overlap = overlap;
661dd8e54a2SToby Isaac   PetscFunctionReturn(0);
662dd8e54a2SToby Isaac }
663dd8e54a2SToby Isaac 
6649be51f97SToby Isaac /*@
6659be51f97SToby Isaac   DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values
6669be51f97SToby Isaac   > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells
6679be51f97SToby Isaac 
6689be51f97SToby Isaac   Not collective
6699be51f97SToby Isaac 
6709be51f97SToby Isaac   Input Parameter:
6719be51f97SToby Isaac . dm - the forest
6729be51f97SToby Isaac 
6739be51f97SToby Isaac   Output Parameter:
6749be51f97SToby Isaac . overlap - default 0
6759be51f97SToby Isaac 
6769be51f97SToby Isaac   Level: intermediate
6779be51f97SToby Isaac 
678db781477SPatrick Sanan .seealso: `DMForestGetPartitionOverlap()`, `DMForestSetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`
6799be51f97SToby Isaac @*/
680dd8e54a2SToby Isaac PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap)
681dd8e54a2SToby Isaac {
682dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
683dd8e54a2SToby Isaac 
684dd8e54a2SToby Isaac   PetscFunctionBegin;
685dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
686dd8e54a2SToby Isaac   PetscValidIntPointer(overlap,2);
687dd8e54a2SToby Isaac   *overlap = forest->overlap;
688dd8e54a2SToby Isaac   PetscFunctionReturn(0);
689dd8e54a2SToby Isaac }
690dd8e54a2SToby Isaac 
6919be51f97SToby Isaac /*@
6929be51f97SToby Isaac   DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base
6939be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest
6949be51f97SToby Isaac   (see DMForestGetAdaptivityForest()) this limits the amount of coarsening.
6959be51f97SToby Isaac 
6969be51f97SToby Isaac   Logically collective on dm
6979be51f97SToby Isaac 
6989be51f97SToby Isaac   Input Parameters:
6999be51f97SToby Isaac + dm - the forest
7009be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7019be51f97SToby Isaac 
7029be51f97SToby Isaac   Level: intermediate
7039be51f97SToby Isaac 
704db781477SPatrick Sanan .seealso: `DMForestGetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestSetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()`
7059be51f97SToby Isaac @*/
706dd8e54a2SToby Isaac PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement)
707dd8e54a2SToby Isaac {
708dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
709dd8e54a2SToby Isaac 
710dd8e54a2SToby Isaac   PetscFunctionBegin;
711dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
71228b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup");
713dd8e54a2SToby Isaac   forest->minRefinement = minRefinement;
714dd8e54a2SToby Isaac   PetscFunctionReturn(0);
715dd8e54a2SToby Isaac }
716dd8e54a2SToby Isaac 
7179be51f97SToby Isaac /*@
7189be51f97SToby Isaac   DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see
7199be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest (see
7209be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of coarsening.
7219be51f97SToby Isaac 
7229be51f97SToby Isaac   Not collective
7239be51f97SToby Isaac 
7249be51f97SToby Isaac   Input Parameter:
7259be51f97SToby Isaac . dm - the forest
7269be51f97SToby Isaac 
7279be51f97SToby Isaac   Output Parameter:
7289be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7299be51f97SToby Isaac 
7309be51f97SToby Isaac   Level: intermediate
7319be51f97SToby Isaac 
732db781477SPatrick Sanan .seealso: `DMForestSetMinimumRefinement()`, `DMForestGetMaximumRefinement()`, `DMForestGetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()`
7339be51f97SToby Isaac @*/
734dd8e54a2SToby Isaac PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement)
735dd8e54a2SToby Isaac {
736dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
737dd8e54a2SToby Isaac 
738dd8e54a2SToby Isaac   PetscFunctionBegin;
739dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
740dd8e54a2SToby Isaac   PetscValidIntPointer(minRefinement,2);
741dd8e54a2SToby Isaac   *minRefinement = forest->minRefinement;
742dd8e54a2SToby Isaac   PetscFunctionReturn(0);
743dd8e54a2SToby Isaac }
744dd8e54a2SToby Isaac 
7459be51f97SToby Isaac /*@
7469be51f97SToby Isaac   DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base
7479be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.
7489be51f97SToby Isaac 
7499be51f97SToby Isaac   Logically collective on dm
7509be51f97SToby Isaac 
7519be51f97SToby Isaac   Input Parameters:
7529be51f97SToby Isaac + dm - the forest
7539be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7549be51f97SToby Isaac 
7559be51f97SToby Isaac   Level: intermediate
7569be51f97SToby Isaac 
757db781477SPatrick Sanan .seealso: `DMForestSetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestGetBaseDM()`
7589be51f97SToby Isaac @*/
75956ba9f64SToby Isaac PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement)
76056ba9f64SToby Isaac {
76156ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
76256ba9f64SToby Isaac 
76356ba9f64SToby Isaac   PetscFunctionBegin;
76456ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
76528b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup");
76656ba9f64SToby Isaac   forest->initRefinement = initRefinement;
76756ba9f64SToby Isaac   PetscFunctionReturn(0);
76856ba9f64SToby Isaac }
76956ba9f64SToby Isaac 
7709be51f97SToby Isaac /*@
7719be51f97SToby Isaac   DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see
7729be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.
7739be51f97SToby Isaac 
7749be51f97SToby Isaac   Not collective
7759be51f97SToby Isaac 
7769be51f97SToby Isaac   Input Parameter:
7779be51f97SToby Isaac . dm - the forest
7789be51f97SToby Isaac 
77901d2d390SJose E. Roman   Output Parameter:
780bf2d5fbbSStefano Zampini . initRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7819be51f97SToby Isaac 
7829be51f97SToby Isaac   Level: intermediate
7839be51f97SToby Isaac 
784db781477SPatrick Sanan .seealso: `DMForestSetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestGetBaseDM()`
7859be51f97SToby Isaac @*/
78656ba9f64SToby Isaac PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement)
78756ba9f64SToby Isaac {
78856ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
78956ba9f64SToby Isaac 
79056ba9f64SToby Isaac   PetscFunctionBegin;
79156ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
79256ba9f64SToby Isaac   PetscValidIntPointer(initRefinement,2);
79356ba9f64SToby Isaac   *initRefinement = forest->initRefinement;
79456ba9f64SToby Isaac   PetscFunctionReturn(0);
79556ba9f64SToby Isaac }
79656ba9f64SToby Isaac 
7979be51f97SToby Isaac /*@
7989be51f97SToby Isaac   DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base
7999be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest
8009be51f97SToby Isaac   (see DMForestGetAdaptivityForest()), this limits the amount of refinement.
8019be51f97SToby Isaac 
8029be51f97SToby Isaac   Logically collective on dm
8039be51f97SToby Isaac 
8049be51f97SToby Isaac   Input Parameters:
8059be51f97SToby Isaac + dm - the forest
8069be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8079be51f97SToby Isaac 
8089be51f97SToby Isaac   Level: intermediate
8099be51f97SToby Isaac 
810db781477SPatrick Sanan .seealso: `DMForestGetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestSetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityDM()`
8119be51f97SToby Isaac @*/
812c7eeac06SToby Isaac PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement)
813dd8e54a2SToby Isaac {
814dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
815dd8e54a2SToby Isaac 
816dd8e54a2SToby Isaac   PetscFunctionBegin;
817dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
81828b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup");
819c7eeac06SToby Isaac   forest->maxRefinement = maxRefinement;
820dd8e54a2SToby Isaac   PetscFunctionReturn(0);
821dd8e54a2SToby Isaac }
822dd8e54a2SToby Isaac 
8239be51f97SToby Isaac /*@
8249be51f97SToby Isaac   DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see
8259be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest (see
8269be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of refinement.
8279be51f97SToby Isaac 
8289be51f97SToby Isaac   Not collective
8299be51f97SToby Isaac 
8309be51f97SToby Isaac   Input Parameter:
8319be51f97SToby Isaac . dm - the forest
8329be51f97SToby Isaac 
8339be51f97SToby Isaac   Output Parameter:
8349be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8359be51f97SToby Isaac 
8369be51f97SToby Isaac   Level: intermediate
8379be51f97SToby Isaac 
838db781477SPatrick Sanan .seealso: `DMForestSetMaximumRefinement()`, `DMForestGetMinimumRefinement()`, `DMForestGetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()`
8399be51f97SToby Isaac @*/
840c7eeac06SToby Isaac PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement)
841dd8e54a2SToby Isaac {
842dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
843dd8e54a2SToby Isaac 
844dd8e54a2SToby Isaac   PetscFunctionBegin;
845dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
846c7eeac06SToby Isaac   PetscValidIntPointer(maxRefinement,2);
847c7eeac06SToby Isaac   *maxRefinement = forest->maxRefinement;
848dd8e54a2SToby Isaac   PetscFunctionReturn(0);
849dd8e54a2SToby Isaac }
850c7eeac06SToby Isaac 
8519be51f97SToby Isaac /*@C
8529be51f97SToby Isaac   DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes.
8539be51f97SToby Isaac   Subtypes of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree
8549be51f97SToby Isaac   for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to
8559be51f97SToby Isaac   specify refinement/coarsening.
8569be51f97SToby Isaac 
8579be51f97SToby Isaac   Logically collective on dm
8589be51f97SToby Isaac 
8599be51f97SToby Isaac   Input Parameters:
8609be51f97SToby Isaac + dm - the forest
8619be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL
8629be51f97SToby Isaac 
8639be51f97SToby Isaac   Level: advanced
8649be51f97SToby Isaac 
865db781477SPatrick Sanan .seealso: `DMForestGetAdaptivityStrategy()`
8669be51f97SToby Isaac @*/
867c7eeac06SToby Isaac PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy)
868c7eeac06SToby Isaac {
869c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
870c7eeac06SToby Isaac 
871c7eeac06SToby Isaac   PetscFunctionBegin;
872c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8739566063dSJacob Faibussowitsch   PetscCall(PetscFree(forest->adaptStrategy));
8749566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy));
875c7eeac06SToby Isaac   PetscFunctionReturn(0);
876c7eeac06SToby Isaac }
877c7eeac06SToby Isaac 
8789be51f97SToby Isaac /*@C
8799be51f97SToby Isaac   DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes.  Subtypes
8809be51f97SToby Isaac   of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all
8819be51f97SToby Isaac   processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only
8829be51f97SToby Isaac   one process needs to specify refinement/coarsening.
8839be51f97SToby Isaac 
8849be51f97SToby Isaac   Not collective
8859be51f97SToby Isaac 
8869be51f97SToby Isaac   Input Parameter:
8879be51f97SToby Isaac . dm - the forest
8889be51f97SToby Isaac 
8899be51f97SToby Isaac   Output Parameter:
8909be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL)
8919be51f97SToby Isaac 
8929be51f97SToby Isaac   Level: advanced
8939be51f97SToby Isaac 
894db781477SPatrick Sanan .seealso: `DMForestSetAdaptivityStrategy()`
8959be51f97SToby Isaac @*/
896c7eeac06SToby Isaac PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
897c7eeac06SToby Isaac {
898c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
899c7eeac06SToby Isaac 
900c7eeac06SToby Isaac   PetscFunctionBegin;
901c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
902c7eeac06SToby Isaac   PetscValidPointer(adaptStrategy,2);
903c7eeac06SToby Isaac   *adaptStrategy = forest->adaptStrategy;
904c7eeac06SToby Isaac   PetscFunctionReturn(0);
905c7eeac06SToby Isaac }
906c7eeac06SToby Isaac 
9072a133e43SToby Isaac /*@
9082a133e43SToby Isaac   DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning,
9092a133e43SToby Isaac   etc.) was successful.  PETSC_FALSE indicates that the post-adaptation forest is the same as the pre-adpatation
9102a133e43SToby Isaac   forest.  A requested adaptation may have been unsuccessful if, for example, the requested refinement would have
9112a133e43SToby Isaac   exceeded the maximum refinement level.
9122a133e43SToby Isaac 
9132a133e43SToby Isaac   Collective on dm
9142a133e43SToby Isaac 
9152a133e43SToby Isaac   Input Parameter:
9162a133e43SToby Isaac 
9172a133e43SToby Isaac . dm - the post-adaptation forest
9182a133e43SToby Isaac 
9192a133e43SToby Isaac   Output Parameter:
9202a133e43SToby Isaac 
9212a133e43SToby Isaac . success - PETSC_TRUE if the post-adaptation forest is different from the pre-adaptation forest.
9222a133e43SToby Isaac 
9232a133e43SToby Isaac   Level: intermediate
9242a133e43SToby Isaac 
9252a133e43SToby Isaac .see
9262a133e43SToby Isaac @*/
9272a133e43SToby Isaac PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success)
9282a133e43SToby Isaac {
9292a133e43SToby Isaac   DM_Forest      *forest;
9302a133e43SToby Isaac 
9312a133e43SToby Isaac   PetscFunctionBegin;
9322a133e43SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
93328b400f6SJacob Faibussowitsch   PetscCheck(dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() has not been called yet.");
9342a133e43SToby Isaac   forest = (DM_Forest *) dm->data;
9359566063dSJacob Faibussowitsch   PetscCall((forest->getadaptivitysuccess)(dm,success));
9362a133e43SToby Isaac   PetscFunctionReturn(0);
9372a133e43SToby Isaac }
9382a133e43SToby Isaac 
939bf9b5d84SToby Isaac /*@
940bf9b5d84SToby Isaac   DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed
941bf9b5d84SToby Isaac   relating the cells of the pre-adaptation forest to the post-adaptiation forest.  After DMSetUp() is called, these transfer PetscSFs can be accessed with DMForestGetAdaptivitySF().
942bf9b5d84SToby Isaac 
943bf9b5d84SToby Isaac   Logically collective on dm
944bf9b5d84SToby Isaac 
945bf9b5d84SToby Isaac   Input Parameters:
946bf9b5d84SToby Isaac + dm - the post-adaptation forest
947bf9b5d84SToby Isaac - computeSF - default PETSC_TRUE
948bf9b5d84SToby Isaac 
949bf9b5d84SToby Isaac   Level: advanced
950bf9b5d84SToby Isaac 
951db781477SPatrick Sanan .seealso: `DMForestGetComputeAdaptivitySF()`, `DMForestGetAdaptivitySF()`
952bf9b5d84SToby Isaac @*/
953bf9b5d84SToby Isaac PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
954bf9b5d84SToby Isaac {
955bf9b5d84SToby Isaac   DM_Forest *forest;
956bf9b5d84SToby Isaac 
957bf9b5d84SToby Isaac   PetscFunctionBegin;
958bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
95928b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called");
960bf9b5d84SToby Isaac   forest                 = (DM_Forest*) dm->data;
961bf9b5d84SToby Isaac   forest->computeAdaptSF = computeSF;
962bf9b5d84SToby Isaac   PetscFunctionReturn(0);
963bf9b5d84SToby Isaac }
964bf9b5d84SToby Isaac 
9650eb7e1eaSToby Isaac PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time)
96680b27e07SToby Isaac {
96780b27e07SToby Isaac   DM_Forest      *forest;
96880b27e07SToby Isaac 
96980b27e07SToby Isaac   PetscFunctionBegin;
97080b27e07SToby Isaac   PetscValidHeaderSpecific(dmIn   ,DM_CLASSID  ,1);
97180b27e07SToby Isaac   PetscValidHeaderSpecific(vecIn  ,VEC_CLASSID ,2);
97280b27e07SToby Isaac   PetscValidHeaderSpecific(dmOut  ,DM_CLASSID  ,3);
97380b27e07SToby Isaac   PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,4);
97480b27e07SToby Isaac   forest = (DM_Forest *) dmIn->data;
97528b400f6SJacob Faibussowitsch   PetscCheck(forest->transfervec,PetscObjectComm((PetscObject)dmIn),PETSC_ERR_SUP,"DMForestTransferVec() not implemented");
9769566063dSJacob Faibussowitsch   PetscCall((forest->transfervec)(dmIn,vecIn,dmOut,vecOut,useBCs,time));
97780b27e07SToby Isaac   PetscFunctionReturn(0);
97880b27e07SToby Isaac }
97980b27e07SToby Isaac 
980ac34a06fSStefano Zampini PetscErrorCode DMForestTransferVecFromBase(DM dm, Vec vecIn, Vec vecOut)
981ac34a06fSStefano Zampini {
982ac34a06fSStefano Zampini   DM_Forest      *forest;
983ac34a06fSStefano Zampini 
984ac34a06fSStefano Zampini   PetscFunctionBegin;
985ac34a06fSStefano Zampini   PetscValidHeaderSpecific(dm   ,DM_CLASSID  ,1);
986ac34a06fSStefano Zampini   PetscValidHeaderSpecific(vecIn  ,VEC_CLASSID ,2);
987ac34a06fSStefano Zampini   PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,3);
988ac34a06fSStefano Zampini   forest = (DM_Forest *) dm->data;
98928b400f6SJacob Faibussowitsch   PetscCheck(forest->transfervecfrombase,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DMForestTransferVecFromBase() not implemented");
9909566063dSJacob Faibussowitsch   PetscCall((forest->transfervecfrombase)(dm,vecIn,vecOut));
991ac34a06fSStefano Zampini   PetscFunctionReturn(0);
992ac34a06fSStefano Zampini }
993ac34a06fSStefano Zampini 
994bf9b5d84SToby Isaac /*@
995bf9b5d84SToby Isaac   DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the
996bf9b5d84SToby Isaac   pre-adaptation forest to the post-adaptiation forest.  After DMSetUp() is called, these transfer PetscSFs can be
997bf9b5d84SToby Isaac   accessed with DMForestGetAdaptivitySF().
998bf9b5d84SToby Isaac 
999bf9b5d84SToby Isaac   Not collective
1000bf9b5d84SToby Isaac 
1001bf9b5d84SToby Isaac   Input Parameter:
1002bf9b5d84SToby Isaac . dm - the post-adaptation forest
1003bf9b5d84SToby Isaac 
1004bf9b5d84SToby Isaac   Output Parameter:
1005bf9b5d84SToby Isaac . computeSF - default PETSC_TRUE
1006bf9b5d84SToby Isaac 
1007bf9b5d84SToby Isaac   Level: advanced
1008bf9b5d84SToby Isaac 
1009db781477SPatrick Sanan .seealso: `DMForestSetComputeAdaptivitySF()`, `DMForestGetAdaptivitySF()`
1010bf9b5d84SToby Isaac @*/
1011bf9b5d84SToby Isaac PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1012bf9b5d84SToby Isaac {
1013bf9b5d84SToby Isaac   DM_Forest *forest;
1014bf9b5d84SToby Isaac 
1015bf9b5d84SToby Isaac   PetscFunctionBegin;
1016bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1017bf9b5d84SToby Isaac   forest     = (DM_Forest*) dm->data;
1018bf9b5d84SToby Isaac   *computeSF = forest->computeAdaptSF;
1019bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1020bf9b5d84SToby Isaac }
1021bf9b5d84SToby Isaac 
1022bf9b5d84SToby Isaac /*@
1023bf9b5d84SToby Isaac   DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest.
1024bf9b5d84SToby Isaac   Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be
1025bf9b5d84SToby Isaac   some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa.  Therefore there are two
1026bf9b5d84SToby Isaac   PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates
1027bf9b5d84SToby Isaac   pre-adaptation fine cells to post-adaptation coarse cells.
1028bf9b5d84SToby Isaac 
1029bf9b5d84SToby Isaac   Not collective
1030bf9b5d84SToby Isaac 
1031bf9b5d84SToby Isaac   Input Parameter:
1032bf9b5d84SToby Isaac   dm - the post-adaptation forest
1033bf9b5d84SToby Isaac 
1034bf9b5d84SToby Isaac   Output Parameter:
10350f17b9e3SToby Isaac   preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
10360f17b9e3SToby Isaac   coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-
1037bf9b5d84SToby Isaac 
1038bf9b5d84SToby Isaac   Level: advanced
1039bf9b5d84SToby Isaac 
1040db781477SPatrick Sanan .seealso: `DMForestGetComputeAdaptivitySF()`, `DMForestSetComputeAdaptivitySF()`
1041bf9b5d84SToby Isaac @*/
10420f17b9e3SToby Isaac PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1043bf9b5d84SToby Isaac {
1044bf9b5d84SToby Isaac   DM_Forest      *forest;
1045bf9b5d84SToby Isaac 
1046bf9b5d84SToby Isaac   PetscFunctionBegin;
1047bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10489566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
1049bf9b5d84SToby Isaac   forest = (DM_Forest*) dm->data;
1050f885a11aSToby Isaac   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1051f885a11aSToby Isaac   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
1052bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1053bf9b5d84SToby Isaac }
1054bf9b5d84SToby Isaac 
10559be51f97SToby Isaac /*@
10569be51f97SToby Isaac   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
10579be51f97SToby Isaac   indicate that the diameter of neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may
10589be51f97SToby Isaac   only support one particular choice of grading factor.
10599be51f97SToby Isaac 
10609be51f97SToby Isaac   Logically collective on dm
10619be51f97SToby Isaac 
10629be51f97SToby Isaac   Input Parameters:
10639be51f97SToby Isaac + dm - the forest
10649be51f97SToby Isaac - grade - the grading factor
10659be51f97SToby Isaac 
10669be51f97SToby Isaac   Level: advanced
10679be51f97SToby Isaac 
1068db781477SPatrick Sanan .seealso: `DMForestGetGradeFactor()`
10699be51f97SToby Isaac @*/
1070c7eeac06SToby Isaac PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1071c7eeac06SToby Isaac {
1072c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1073c7eeac06SToby Isaac 
1074c7eeac06SToby Isaac   PetscFunctionBegin;
1075c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
107628b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup");
1077c7eeac06SToby Isaac   forest->gradeFactor = grade;
1078c7eeac06SToby Isaac   PetscFunctionReturn(0);
1079c7eeac06SToby Isaac }
1080c7eeac06SToby Isaac 
10819be51f97SToby Isaac /*@
10829be51f97SToby Isaac   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
10839be51f97SToby Isaac   neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may only support one particular
10849be51f97SToby Isaac   choice of grading factor.
10859be51f97SToby Isaac 
10869be51f97SToby Isaac   Not collective
10879be51f97SToby Isaac 
10889be51f97SToby Isaac   Input Parameter:
10899be51f97SToby Isaac . dm - the forest
10909be51f97SToby Isaac 
10919be51f97SToby Isaac   Output Parameter:
10929be51f97SToby Isaac . grade - the grading factor
10939be51f97SToby Isaac 
10949be51f97SToby Isaac   Level: advanced
10959be51f97SToby Isaac 
1096db781477SPatrick Sanan .seealso: `DMForestSetGradeFactor()`
10979be51f97SToby Isaac @*/
1098c7eeac06SToby Isaac PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1099c7eeac06SToby Isaac {
1100c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1101c7eeac06SToby Isaac 
1102c7eeac06SToby Isaac   PetscFunctionBegin;
1103c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1104c7eeac06SToby Isaac   PetscValidIntPointer(grade,2);
1105c7eeac06SToby Isaac   *grade = forest->gradeFactor;
1106c7eeac06SToby Isaac   PetscFunctionReturn(0);
1107c7eeac06SToby Isaac }
1108c7eeac06SToby Isaac 
11099be51f97SToby Isaac /*@
11109be51f97SToby Isaac   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
11119be51f97SToby Isaac   the cell weight (see DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be
11129be51f97SToby Isaac   (cellWeight) * (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on
11139be51f97SToby Isaac   its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
11149be51f97SToby Isaac   computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11159be51f97SToby Isaac 
11169be51f97SToby Isaac   Logically collective on dm
11179be51f97SToby Isaac 
11189be51f97SToby Isaac   Input Parameters:
11199be51f97SToby Isaac + dm - the forest
11209be51f97SToby Isaac - weightsFactors - default 1.
11219be51f97SToby Isaac 
11229be51f97SToby Isaac   Level: advanced
11239be51f97SToby Isaac 
1124db781477SPatrick Sanan .seealso: `DMForestGetCellWeightFactor()`, `DMForestSetCellWeights()`
11259be51f97SToby Isaac @*/
1126ef51cf95SToby Isaac PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
1127c7eeac06SToby Isaac {
1128c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1129c7eeac06SToby Isaac 
1130c7eeac06SToby Isaac   PetscFunctionBegin;
1131c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
113228b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup");
1133c7eeac06SToby Isaac   forest->weightsFactor = weightsFactor;
1134c7eeac06SToby Isaac   PetscFunctionReturn(0);
1135c7eeac06SToby Isaac }
1136c7eeac06SToby Isaac 
11379be51f97SToby Isaac /*@
11389be51f97SToby Isaac   DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
11399be51f97SToby Isaac   DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be (cellWeight) *
11409be51f97SToby Isaac   (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on its level; a
11419be51f97SToby Isaac   factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation
11429be51f97SToby Isaac   associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11439be51f97SToby Isaac 
11449be51f97SToby Isaac   Not collective
11459be51f97SToby Isaac 
11469be51f97SToby Isaac   Input Parameter:
11479be51f97SToby Isaac . dm - the forest
11489be51f97SToby Isaac 
11499be51f97SToby Isaac   Output Parameter:
11509be51f97SToby Isaac . weightsFactors - default 1.
11519be51f97SToby Isaac 
11529be51f97SToby Isaac   Level: advanced
11539be51f97SToby Isaac 
1154db781477SPatrick Sanan .seealso: `DMForestSetCellWeightFactor()`, `DMForestSetCellWeights()`
11559be51f97SToby Isaac @*/
1156ef51cf95SToby Isaac PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1157c7eeac06SToby Isaac {
1158c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1159c7eeac06SToby Isaac 
1160c7eeac06SToby Isaac   PetscFunctionBegin;
1161c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1162c7eeac06SToby Isaac   PetscValidRealPointer(weightsFactor,2);
1163c7eeac06SToby Isaac   *weightsFactor = forest->weightsFactor;
1164c7eeac06SToby Isaac   PetscFunctionReturn(0);
1165c7eeac06SToby Isaac }
1166c7eeac06SToby Isaac 
11679be51f97SToby Isaac /*@
11689be51f97SToby Isaac   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
11699be51f97SToby Isaac 
11709be51f97SToby Isaac   Not collective
11719be51f97SToby Isaac 
11729be51f97SToby Isaac   Input Parameter:
11739be51f97SToby Isaac . dm - the forest
11749be51f97SToby Isaac 
11759be51f97SToby Isaac   Output Parameters:
11769be51f97SToby Isaac + cStart - the first cell on this process
11779be51f97SToby Isaac - cEnd - one after the final cell on this process
11789be51f97SToby Isaac 
11791a244344SSatish Balay   Level: intermediate
11809be51f97SToby Isaac 
1181db781477SPatrick Sanan .seealso: `DMForestGetCellSF()`
11829be51f97SToby Isaac @*/
1183c7eeac06SToby Isaac PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1184c7eeac06SToby Isaac {
1185c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1186c7eeac06SToby Isaac 
1187c7eeac06SToby Isaac   PetscFunctionBegin;
1188c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1189c7eeac06SToby Isaac   PetscValidIntPointer(cStart,2);
1190064a246eSJacob Faibussowitsch   PetscValidIntPointer(cEnd,3);
1191c7eeac06SToby Isaac   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) {
11929566063dSJacob Faibussowitsch     PetscCall(forest->createcellchart(dm,&forest->cStart,&forest->cEnd));
1193c7eeac06SToby Isaac   }
1194c7eeac06SToby Isaac   *cStart =  forest->cStart;
1195c7eeac06SToby Isaac   *cEnd   =  forest->cEnd;
1196c7eeac06SToby Isaac   PetscFunctionReturn(0);
1197c7eeac06SToby Isaac }
1198c7eeac06SToby Isaac 
11999be51f97SToby Isaac /*@
12009be51f97SToby Isaac   DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes
12019be51f97SToby Isaac 
12029be51f97SToby Isaac   Not collective
12039be51f97SToby Isaac 
12049be51f97SToby Isaac   Input Parameter:
12059be51f97SToby Isaac . dm - the forest
12069be51f97SToby Isaac 
12079be51f97SToby Isaac   Output Parameter:
12089be51f97SToby Isaac . cellSF - the PetscSF
12099be51f97SToby Isaac 
12101a244344SSatish Balay   Level: intermediate
12119be51f97SToby Isaac 
1212db781477SPatrick Sanan .seealso: `DMForestGetCellChart()`
12139be51f97SToby Isaac @*/
1214c7eeac06SToby Isaac PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1215c7eeac06SToby Isaac {
1216c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1217c7eeac06SToby Isaac 
1218c7eeac06SToby Isaac   PetscFunctionBegin;
1219c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1220c7eeac06SToby Isaac   PetscValidPointer(cellSF,2);
1221c7eeac06SToby Isaac   if ((!forest->cellSF) && forest->createcellsf) {
12229566063dSJacob Faibussowitsch     PetscCall(forest->createcellsf(dm,&forest->cellSF));
1223c7eeac06SToby Isaac   }
1224c7eeac06SToby Isaac   *cellSF = forest->cellSF;
1225c7eeac06SToby Isaac   PetscFunctionReturn(0);
1226c7eeac06SToby Isaac }
1227c7eeac06SToby Isaac 
12289be51f97SToby Isaac /*@C
12299be51f97SToby Isaac   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
12309be51f97SToby Isaac   DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination).  The
1231cd3c525cSToby Isaac   interpretation of the label values is up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP,
1232cd3c525cSToby Isaac   DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have been reserved as choices that should be accepted by all subtypes.
12339be51f97SToby Isaac 
12349be51f97SToby Isaac   Logically collective on dm
12359be51f97SToby Isaac 
12369be51f97SToby Isaac   Input Parameters:
12379be51f97SToby Isaac - dm - the forest
1238a1b0c543SToby Isaac + adaptLabel - the label in the pre-adaptation forest
12399be51f97SToby Isaac 
12409be51f97SToby Isaac   Level: intermediate
12419be51f97SToby Isaac 
1242db781477SPatrick Sanan .seealso `DMForestGetAdaptivityLabel()`
12439be51f97SToby Isaac @*/
1244a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel)
1245c7eeac06SToby Isaac {
1246c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1247c7eeac06SToby Isaac 
1248c7eeac06SToby Isaac   PetscFunctionBegin;
1249c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12501fd50544SStefano Zampini   if (adaptLabel) PetscValidHeaderSpecific(adaptLabel,DMLABEL_CLASSID,2);
12519566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)adaptLabel));
12529566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&forest->adaptLabel));
12531fd50544SStefano Zampini   forest->adaptLabel = adaptLabel;
1254c7eeac06SToby Isaac   PetscFunctionReturn(0);
1255c7eeac06SToby Isaac }
1256c7eeac06SToby Isaac 
12579be51f97SToby Isaac /*@C
12589be51f97SToby Isaac   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that
12599be51f97SToby Isaac   holds the adaptation flags (refinement, coarsening, or some combination).  The interpretation of the label values is
1260cd3c525cSToby Isaac   up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have
1261cd3c525cSToby Isaac   been reserved as choices that should be accepted by all subtypes.
12629be51f97SToby Isaac 
12639be51f97SToby Isaac   Not collective
12649be51f97SToby Isaac 
12659be51f97SToby Isaac   Input Parameter:
12669be51f97SToby Isaac . dm - the forest
12679be51f97SToby Isaac 
12689be51f97SToby Isaac   Output Parameter:
12699be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest
12709be51f97SToby Isaac 
12719be51f97SToby Isaac   Level: intermediate
12729be51f97SToby Isaac 
1273db781477SPatrick Sanan .seealso `DMForestSetAdaptivityLabel()`
12749be51f97SToby Isaac @*/
1275a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel)
1276c7eeac06SToby Isaac {
1277c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1278c7eeac06SToby Isaac 
1279c7eeac06SToby Isaac   PetscFunctionBegin;
1280c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1281ba936b91SToby Isaac   *adaptLabel = forest->adaptLabel;
1282c7eeac06SToby Isaac   PetscFunctionReturn(0);
1283c7eeac06SToby Isaac }
1284c7eeac06SToby Isaac 
12859be51f97SToby Isaac /*@
12869be51f97SToby Isaac   DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
12879be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
12889be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
12899be51f97SToby Isaac   of 1.
12909be51f97SToby Isaac 
12919be51f97SToby Isaac   Logically collective on dm
12929be51f97SToby Isaac 
12939be51f97SToby Isaac   Input Parameters:
12949be51f97SToby Isaac + dm - the forest
12959be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
12969be51f97SToby Isaac - copyMode - how weights should reference weights
12979be51f97SToby Isaac 
12989be51f97SToby Isaac   Level: advanced
12999be51f97SToby Isaac 
1300db781477SPatrick Sanan .seealso: `DMForestGetCellWeights()`, `DMForestSetWeightCapacity()`
13019be51f97SToby Isaac @*/
1302c7eeac06SToby Isaac PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1303c7eeac06SToby Isaac {
1304c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1305c7eeac06SToby Isaac   PetscInt       cStart, cEnd;
1306c7eeac06SToby Isaac 
1307c7eeac06SToby Isaac   PetscFunctionBegin;
1308c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
13099566063dSJacob Faibussowitsch   PetscCall(DMForestGetCellChart(dm,&cStart,&cEnd));
131063a3b9bcSJacob Faibussowitsch   PetscCheck(cEnd >= cStart,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%" PetscInt_FMT ",%" PetscInt_FMT ") is not valid",cStart,cEnd);
1311c7eeac06SToby Isaac   if (copyMode == PETSC_COPY_VALUES) {
1312c7eeac06SToby Isaac     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) {
13139566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(cEnd-cStart,&forest->cellWeights));
1314c7eeac06SToby Isaac     }
13159566063dSJacob Faibussowitsch     PetscCall(PetscArraycpy(forest->cellWeights,weights,cEnd-cStart));
1316c7eeac06SToby Isaac     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
1317c7eeac06SToby Isaac     PetscFunctionReturn(0);
1318c7eeac06SToby Isaac   }
1319c7eeac06SToby Isaac   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) {
13209566063dSJacob Faibussowitsch     PetscCall(PetscFree(forest->cellWeights));
1321c7eeac06SToby Isaac   }
1322c7eeac06SToby Isaac   forest->cellWeights         = weights;
1323c7eeac06SToby Isaac   forest->cellWeightsCopyMode = copyMode;
1324c7eeac06SToby Isaac   PetscFunctionReturn(0);
1325c7eeac06SToby Isaac }
1326c7eeac06SToby Isaac 
13279be51f97SToby Isaac /*@
13289be51f97SToby Isaac   DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
13299be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
13309be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
13319be51f97SToby Isaac   of 1.
13329be51f97SToby Isaac 
13339be51f97SToby Isaac   Not collective
13349be51f97SToby Isaac 
13359be51f97SToby Isaac   Input Parameter:
13369be51f97SToby Isaac . dm - the forest
13379be51f97SToby Isaac 
13389be51f97SToby Isaac   Output Parameter:
13399be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
13409be51f97SToby Isaac 
13419be51f97SToby Isaac   Level: advanced
13429be51f97SToby Isaac 
1343db781477SPatrick Sanan .seealso: `DMForestSetCellWeights()`, `DMForestSetWeightCapacity()`
13449be51f97SToby Isaac @*/
1345c7eeac06SToby Isaac PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1346c7eeac06SToby Isaac {
1347c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1348c7eeac06SToby Isaac 
1349c7eeac06SToby Isaac   PetscFunctionBegin;
1350c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1351c7eeac06SToby Isaac   PetscValidPointer(weights,2);
1352c7eeac06SToby Isaac   *weights = forest->cellWeights;
1353c7eeac06SToby Isaac   PetscFunctionReturn(0);
1354c7eeac06SToby Isaac }
1355c7eeac06SToby Isaac 
13569be51f97SToby Isaac /*@
13579be51f97SToby Isaac   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
13589be51f97SToby Isaac   a pre-adaptation forest (see DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each
13599be51f97SToby Isaac   process's cells to the process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that
13609be51f97SToby Isaac   the current process should not have any cells after repartitioning.
13619be51f97SToby Isaac 
13629be51f97SToby Isaac   Logically Collective on dm
13639be51f97SToby Isaac 
13649be51f97SToby Isaac   Input parameters:
13659be51f97SToby Isaac + dm - the forest
13669be51f97SToby Isaac - capacity - this process's capacity
13679be51f97SToby Isaac 
13689be51f97SToby Isaac   Level: advanced
13699be51f97SToby Isaac 
1370db781477SPatrick Sanan .seealso `DMForestGetWeightCapacity()`, `DMForestSetCellWeights()`, `DMForestSetCellWeightFactor()`
13719be51f97SToby Isaac @*/
1372c7eeac06SToby Isaac PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1373c7eeac06SToby Isaac {
1374c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1375c7eeac06SToby Isaac 
1376c7eeac06SToby Isaac   PetscFunctionBegin;
1377c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
137828b400f6SJacob Faibussowitsch   PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup");
137963a3b9bcSJacob Faibussowitsch   PetscCheck(capacity >= 0.,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %g",(double)capacity);
1380c7eeac06SToby Isaac   forest->weightCapacity = capacity;
1381c7eeac06SToby Isaac   PetscFunctionReturn(0);
1382c7eeac06SToby Isaac }
1383c7eeac06SToby Isaac 
13849be51f97SToby Isaac /*@
13859be51f97SToby Isaac   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
13869be51f97SToby Isaac   DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each process's cells to the
13879be51f97SToby Isaac   process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that the current process
13889be51f97SToby Isaac   should not have any cells after repartitioning.
13899be51f97SToby Isaac 
13909be51f97SToby Isaac   Not collective
13919be51f97SToby Isaac 
13929be51f97SToby Isaac   Input parameter:
13939be51f97SToby Isaac . dm - the forest
13949be51f97SToby Isaac 
13959be51f97SToby Isaac   Output parameter:
13969be51f97SToby Isaac . capacity - this process's capacity
13979be51f97SToby Isaac 
13989be51f97SToby Isaac   Level: advanced
13999be51f97SToby Isaac 
1400db781477SPatrick Sanan .seealso `DMForestSetWeightCapacity()`, `DMForestSetCellWeights()`, `DMForestSetCellWeightFactor()`
14019be51f97SToby Isaac @*/
1402c7eeac06SToby Isaac PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1403c7eeac06SToby Isaac {
1404c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1405c7eeac06SToby Isaac 
1406c7eeac06SToby Isaac   PetscFunctionBegin;
1407c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1408c7eeac06SToby Isaac   PetscValidRealPointer(capacity,2);
1409c7eeac06SToby Isaac   *capacity = forest->weightCapacity;
1410c7eeac06SToby Isaac   PetscFunctionReturn(0);
1411c7eeac06SToby Isaac }
1412c7eeac06SToby Isaac 
14130709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm)
1414db4d5e8cSToby Isaac {
141556ba9f64SToby Isaac   PetscBool                  flg, flg1, flg2, flg3, flg4;
1416dd8e54a2SToby Isaac   DMForestTopology           oldTopo;
1417c7eeac06SToby Isaac   char                       stringBuffer[256];
1418dd8e54a2SToby Isaac   PetscViewer                viewer;
1419dd8e54a2SToby Isaac   PetscViewerFormat          format;
142056ba9f64SToby Isaac   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1421c7eeac06SToby Isaac   PetscReal                  weightsFactor;
1422c7eeac06SToby Isaac   DMForestAdaptivityStrategy adaptStrategy;
1423db4d5e8cSToby Isaac 
1424db4d5e8cSToby Isaac   PetscFunctionBegin;
1425db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14269566063dSJacob Faibussowitsch   PetscCall(DMForestGetTopology(dm, &oldTopo));
1427d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"DMForest Options");
14289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,sizeof(stringBuffer),&flg1));
14299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2));
14309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3));
14319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4));
14321dca8a05SBarry 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}");
143356ba9f64SToby Isaac   if (flg1) {
14349566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm,(DMForestTopology)stringBuffer));
14359566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm,NULL));
14369566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm,NULL));
143756ba9f64SToby Isaac   }
143856ba9f64SToby Isaac   if (flg2) {
1439dd8e54a2SToby Isaac     DM base;
1440dd8e54a2SToby Isaac 
14419566063dSJacob Faibussowitsch     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm),&base));
14429566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer,format));
14439566063dSJacob Faibussowitsch     PetscCall(DMLoad(base,viewer));
14449566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
14459566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm,base));
14469566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&base));
14479566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm,NULL));
14489566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm,NULL));
1449dd8e54a2SToby Isaac   }
145056ba9f64SToby Isaac   if (flg3) {
1451dd8e54a2SToby Isaac     DM coarse;
1452dd8e54a2SToby Isaac 
14539566063dSJacob Faibussowitsch     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm),&coarse));
14549566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer,format));
14559566063dSJacob Faibussowitsch     PetscCall(DMLoad(coarse,viewer));
14569566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
14579566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm,coarse));
14589566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&coarse));
14599566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm,NULL));
14609566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm,NULL));
1461dd8e54a2SToby Isaac   }
146256ba9f64SToby Isaac   if (flg4) {
1463dd8e54a2SToby Isaac     DM fine;
1464dd8e54a2SToby Isaac 
14659566063dSJacob Faibussowitsch     PetscCall(DMCreate(PetscObjectComm((PetscObject)dm),&fine));
14669566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer,format));
14679566063dSJacob Faibussowitsch     PetscCall(DMLoad(fine,viewer));
14689566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
14699566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdaptivityForest(dm,fine));
14709566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&fine));
14719566063dSJacob Faibussowitsch     PetscCall(DMForestSetTopology(dm,NULL));
14729566063dSJacob Faibussowitsch     PetscCall(DMForestSetBaseDM(dm,NULL));
1473dd8e54a2SToby Isaac   }
14749566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdjacencyDimension(dm,&adjDim));
14759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg,0));
1476dd8e54a2SToby Isaac   if (flg) {
14779566063dSJacob Faibussowitsch     PetscCall(DMForestSetAdjacencyDimension(dm,adjDim));
1478f885a11aSToby Isaac   } else {
14799566063dSJacob Faibussowitsch     PetscCall(DMForestGetAdjacencyCodimension(dm,&adjCodim));
14809566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg,1));
14811baa6e33SBarry Smith     if (flg) PetscCall(DMForestSetAdjacencyCodimension(dm,adjCodim));
1482dd8e54a2SToby Isaac   }
14839566063dSJacob Faibussowitsch   PetscCall(DMForestGetPartitionOverlap(dm,&overlap));
14849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg,0));
14851baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetPartitionOverlap(dm,overlap));
1486a6121fbdSMatthew G. Knepley #if 0
14879566063dSJacob 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));
1488a6121fbdSMatthew G. Knepley   if (flg) {
14899566063dSJacob Faibussowitsch     PetscCall(DMForestSetMinimumRefinement(dm,minRefinement));
14909566063dSJacob Faibussowitsch     PetscCall(DMForestSetInitialRefinement(dm,minRefinement));
1491a6121fbdSMatthew G. Knepley   }
14929566063dSJacob 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));
1493a6121fbdSMatthew G. Knepley   if (flg) {
14949566063dSJacob Faibussowitsch     PetscCall(DMForestSetMinimumRefinement(dm,0));
14959566063dSJacob Faibussowitsch     PetscCall(DMForestSetInitialRefinement(dm,initRefinement));
1496a6121fbdSMatthew G. Knepley   }
1497a6121fbdSMatthew G. Knepley #endif
14989566063dSJacob Faibussowitsch   PetscCall(DMForestGetMinimumRefinement(dm,&minRefinement));
14999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg,0));
15001baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetMinimumRefinement(dm,minRefinement));
15019566063dSJacob Faibussowitsch   PetscCall(DMForestGetInitialRefinement(dm,&initRefinement));
15029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg,0));
15031baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetInitialRefinement(dm,initRefinement));
15049566063dSJacob Faibussowitsch   PetscCall(DMForestGetMaximumRefinement(dm,&maxRefinement));
15059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg,0));
15061baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetMaximumRefinement(dm,maxRefinement));
15079566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivityStrategy(dm,&adaptStrategy));
15089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,sizeof(stringBuffer),&flg));
15091baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer));
15109566063dSJacob Faibussowitsch   PetscCall(DMForestGetGradeFactor(dm,&grade));
15119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg,0));
15121baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetGradeFactor(dm,grade));
15139566063dSJacob Faibussowitsch   PetscCall(DMForestGetCellWeightFactor(dm,&weightsFactor));
15149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg));
15151baa6e33SBarry Smith   if (flg) PetscCall(DMForestSetCellWeightFactor(dm,weightsFactor));
1516d0609cedSBarry Smith   PetscOptionsHeadEnd();
1517db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1518db4d5e8cSToby Isaac }
1519db4d5e8cSToby Isaac 
1520276c5506SMatthew G. Knepley PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1521d8984e3bSMatthew G. Knepley {
1522d8984e3bSMatthew G. Knepley   PetscFunctionBegin;
15239566063dSJacob Faibussowitsch   if (subdm) PetscCall(DMClone(dm, subdm));
15249566063dSJacob Faibussowitsch   PetscCall(DMCreateSectionSubDM(dm, numFields, fields, is, subdm));
1525d8984e3bSMatthew G. Knepley   PetscFunctionReturn(0);
1526d8984e3bSMatthew G. Knepley }
1527d8984e3bSMatthew G. Knepley 
15285421bac9SToby Isaac PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
15295421bac9SToby Isaac {
15305421bac9SToby Isaac   DMLabel        refine;
15315421bac9SToby Isaac   DM             fineDM;
15325421bac9SToby Isaac 
15335421bac9SToby Isaac   PetscFunctionBegin;
15349566063dSJacob Faibussowitsch   PetscCall(DMGetFineDM(dm,&fineDM));
15355421bac9SToby Isaac   if (fineDM) {
15369566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)fineDM));
15375421bac9SToby Isaac     *dmRefined = fineDM;
15385421bac9SToby Isaac     PetscFunctionReturn(0);
15395421bac9SToby Isaac   }
15409566063dSJacob Faibussowitsch   PetscCall(DMForestTemplate(dm,comm,dmRefined));
15419566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm,"refine",&refine));
15425421bac9SToby Isaac   if (!refine) {
15439566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "refine",&refine));
15449566063dSJacob Faibussowitsch     PetscCall(DMLabelSetDefaultValue(refine,DM_ADAPT_REFINE));
15451baa6e33SBarry Smith   } else PetscCall(PetscObjectReference((PetscObject) refine));
15469566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityLabel(*dmRefined,refine));
15479566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&refine));
15485421bac9SToby Isaac   PetscFunctionReturn(0);
15495421bac9SToby Isaac }
15505421bac9SToby Isaac 
15515421bac9SToby Isaac PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
15525421bac9SToby Isaac {
15535421bac9SToby Isaac   DMLabel        coarsen;
15545421bac9SToby Isaac   DM             coarseDM;
15555421bac9SToby Isaac 
15565421bac9SToby Isaac   PetscFunctionBegin;
15574098eed7SToby Isaac   {
15584098eed7SToby Isaac     PetscMPIInt mpiComparison;
15594098eed7SToby Isaac     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);
15604098eed7SToby Isaac 
15619566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_compare(comm, dmcomm, &mpiComparison));
15621dca8a05SBarry Smith     PetscCheck(mpiComparison == MPI_IDENT || mpiComparison == MPI_CONGRUENT,dmcomm,PETSC_ERR_SUP,"No support for different communicators yet");
15634098eed7SToby Isaac   }
15649566063dSJacob Faibussowitsch   PetscCall(DMGetCoarseDM(dm,&coarseDM));
15655421bac9SToby Isaac   if (coarseDM) {
15669566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)coarseDM));
15675421bac9SToby Isaac     *dmCoarsened = coarseDM;
15685421bac9SToby Isaac     PetscFunctionReturn(0);
15695421bac9SToby Isaac   }
15709566063dSJacob Faibussowitsch   PetscCall(DMForestTemplate(dm,comm,dmCoarsened));
15719566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityPurpose(*dmCoarsened,DM_ADAPT_COARSEN));
15729566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm,"coarsen",&coarsen));
15735421bac9SToby Isaac   if (!coarsen) {
15749566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "coarsen",&coarsen));
15759566063dSJacob Faibussowitsch     PetscCall(DMLabelSetDefaultValue(coarsen,DM_ADAPT_COARSEN));
15761baa6e33SBarry Smith   } else PetscCall(PetscObjectReference((PetscObject) coarsen));
15779566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityLabel(*dmCoarsened,coarsen));
15789566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&coarsen));
15795421bac9SToby Isaac   PetscFunctionReturn(0);
15805421bac9SToby Isaac }
15815421bac9SToby Isaac 
15829fe9e680SJoe Wallwork PetscErrorCode DMAdaptLabel_Forest(DM dm, PETSC_UNUSED Vec metric, DMLabel label, PETSC_UNUSED DMLabel rgLabel, DM *adaptedDM)
158309350103SToby Isaac {
158409350103SToby Isaac   PetscBool      success;
158509350103SToby Isaac 
158609350103SToby Isaac   PetscFunctionBegin;
15879566063dSJacob Faibussowitsch   PetscCall(DMForestTemplate(dm,PetscObjectComm((PetscObject)dm),adaptedDM));
15889566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityLabel(*adaptedDM,label));
15899566063dSJacob Faibussowitsch   PetscCall(DMSetUp(*adaptedDM));
15909566063dSJacob Faibussowitsch   PetscCall(DMForestGetAdaptivitySuccess(*adaptedDM,&success));
159109350103SToby Isaac   if (!success) {
15929566063dSJacob Faibussowitsch     PetscCall(DMDestroy(adaptedDM));
159309350103SToby Isaac     *adaptedDM = NULL;
159409350103SToby Isaac   }
159509350103SToby Isaac   PetscFunctionReturn(0);
159609350103SToby Isaac }
159709350103SToby Isaac 
1598d222f98bSToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm)
1599d222f98bSToby Isaac {
1600d222f98bSToby Isaac   PetscFunctionBegin;
16019566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(dm->ops,sizeof(*(dm->ops))));
1602d222f98bSToby Isaac 
1603d222f98bSToby Isaac   dm->ops->clone          = DMClone_Forest;
1604d222f98bSToby Isaac   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1605d222f98bSToby Isaac   dm->ops->destroy        = DMDestroy_Forest;
1606d8984e3bSMatthew G. Knepley   dm->ops->createsubdm    = DMCreateSubDM_Forest;
16075421bac9SToby Isaac   dm->ops->refine         = DMRefine_Forest;
16085421bac9SToby Isaac   dm->ops->coarsen        = DMCoarsen_Forest;
1609d222f98bSToby Isaac   PetscFunctionReturn(0);
1610d222f98bSToby Isaac }
1611d222f98bSToby Isaac 
16129be51f97SToby Isaac /*MC
16139be51f97SToby Isaac 
1614bae1f979SBarry Smith      DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh.  Forests usually have a base DM
1615bae1f979SBarry Smith   (see DMForestGetBaseDM()), from which it is refined.  The refinement and partitioning of forests is considered
1616bae1f979SBarry Smith   immutable after DMSetUp() is called.  To adapt a mesh, one should call DMForestTemplate() to create a new mesh that
1617bae1f979SBarry Smith   will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new
1618bae1f979SBarry Smith   mesh.
1619bae1f979SBarry Smith 
1620bae1f979SBarry Smith   To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the
1621bae1f979SBarry Smith   previous mesh whose values indicate which cells should be refined (DM_ADAPT_REFINE) or coarsened (DM_ADAPT_COARSEN)
1622bae1f979SBarry Smith   and how (subtypes are free to allow additional values for things like anisotropic refinement).  The label should be
1623bae1f979SBarry Smith   given to the *new* mesh with DMForestSetAdaptivityLabel().
16249be51f97SToby Isaac 
16259be51f97SToby Isaac   Level: advanced
16269be51f97SToby Isaac 
1627db781477SPatrick Sanan .seealso: `DMType`, `DMCreate()`, `DMSetType()`, `DMForestGetBaseDM()`, `DMForestSetBaseDM()`, `DMForestTemplate()`, `DMForestSetAdaptivityLabel()`
16289be51f97SToby Isaac M*/
16299be51f97SToby Isaac 
1630db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1631db4d5e8cSToby Isaac {
1632db4d5e8cSToby Isaac   DM_Forest      *forest;
1633db4d5e8cSToby Isaac 
1634db4d5e8cSToby Isaac   PetscFunctionBegin;
1635db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
16369566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(dm,&forest));
1637db4d5e8cSToby Isaac   dm->dim                      = 0;
1638db4d5e8cSToby Isaac   dm->data                     = forest;
1639db4d5e8cSToby Isaac   forest->refct                = 1;
1640db4d5e8cSToby Isaac   forest->data                 = NULL;
1641db4d5e8cSToby Isaac   forest->topology             = NULL;
1642cd3c525cSToby Isaac   forest->adapt                = NULL;
1643db4d5e8cSToby Isaac   forest->base                 = NULL;
16446a87ffbfSToby Isaac   forest->adaptPurpose         = DM_ADAPT_DETERMINE;
1645db4d5e8cSToby Isaac   forest->adjDim               = PETSC_DEFAULT;
1646db4d5e8cSToby Isaac   forest->overlap              = PETSC_DEFAULT;
1647db4d5e8cSToby Isaac   forest->minRefinement        = PETSC_DEFAULT;
1648db4d5e8cSToby Isaac   forest->maxRefinement        = PETSC_DEFAULT;
164956ba9f64SToby Isaac   forest->initRefinement       = PETSC_DEFAULT;
1650c7eeac06SToby Isaac   forest->cStart               = PETSC_DETERMINE;
1651c7eeac06SToby Isaac   forest->cEnd                 = PETSC_DETERMINE;
1652cd3c525cSToby Isaac   forest->cellSF               = NULL;
1653ebdf65a2SToby Isaac   forest->adaptLabel           = NULL;
1654db4d5e8cSToby Isaac   forest->gradeFactor          = 2;
1655db4d5e8cSToby Isaac   forest->cellWeights          = NULL;
1656db4d5e8cSToby Isaac   forest->cellWeightsCopyMode  = PETSC_USE_POINTER;
1657db4d5e8cSToby Isaac   forest->weightsFactor        = 1.;
1658db4d5e8cSToby Isaac   forest->weightCapacity       = 1.;
16599566063dSJacob Faibussowitsch   PetscCall(DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL));
16609566063dSJacob Faibussowitsch   PetscCall(DMInitialize_Forest(dm));
1661db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1662db4d5e8cSToby Isaac }
1663