xref: /petsc/src/dm/impls/forest/forest.c (revision 456cc5b74ec521bbe64a67e1fc66cb6e7984bbda)
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   PetscErrorCode   ierr;
2227d4645fSToby Isaac 
2327d4645fSToby Isaac   PetscFunctionBegin;
2427d4645fSToby Isaac   while (link) {
2527d4645fSToby Isaac     oldLink = link;
26f416af30SBarry Smith     ierr    = PetscFree(oldLink->name);CHKERRQ(ierr);
2727d4645fSToby Isaac     link    = oldLink->next;
2827d4645fSToby Isaac     ierr    = PetscFree(oldLink);CHKERRQ(ierr);
2927d4645fSToby Isaac   }
3027d4645fSToby Isaac   PetscFunctionReturn(0);
3127d4645fSToby Isaac }
3227d4645fSToby Isaac 
3327d4645fSToby Isaac static PetscErrorCode DMForestPackageInitialize(void)
3427d4645fSToby Isaac {
3527d4645fSToby Isaac   PetscErrorCode ierr;
3627d4645fSToby Isaac 
3727d4645fSToby Isaac   PetscFunctionBegin;
3827d4645fSToby Isaac   if (DMForestPackageInitialized) PetscFunctionReturn(0);
3927d4645fSToby Isaac   DMForestPackageInitialized = PETSC_TRUE;
40f885a11aSToby Isaac 
4127d4645fSToby Isaac   ierr = DMForestRegisterType(DMFOREST);CHKERRQ(ierr);
4227d4645fSToby Isaac   ierr = PetscRegisterFinalize(DMForestPackageFinalize);CHKERRQ(ierr);
4327d4645fSToby Isaac   PetscFunctionReturn(0);
4427d4645fSToby Isaac }
4527d4645fSToby Isaac 
469be51f97SToby Isaac /*@C
479be51f97SToby Isaac   DMForestRegisterType - Registers a DMType as a subtype of DMFOREST (so that DMIsForest() will be correct)
489be51f97SToby Isaac 
499be51f97SToby Isaac   Not Collective
509be51f97SToby Isaac 
519be51f97SToby Isaac   Input parameter:
529be51f97SToby Isaac . name - the name of the type
539be51f97SToby Isaac 
549be51f97SToby Isaac   Level: advanced
559be51f97SToby Isaac 
569be51f97SToby Isaac .seealso: DMFOREST, DMIsForest()
579be51f97SToby Isaac @*/
5827d4645fSToby Isaac PetscErrorCode DMForestRegisterType(DMType name)
5927d4645fSToby Isaac {
6027d4645fSToby Isaac   DMForestTypeLink link;
6127d4645fSToby Isaac   PetscErrorCode   ierr;
6227d4645fSToby Isaac 
6327d4645fSToby Isaac   PetscFunctionBegin;
6427d4645fSToby Isaac   ierr             = DMForestPackageInitialize();CHKERRQ(ierr);
6527d4645fSToby Isaac   ierr             = PetscNew(&link);CHKERRQ(ierr);
6627d4645fSToby Isaac   ierr             = PetscStrallocpy(name,&link->name);CHKERRQ(ierr);
6727d4645fSToby Isaac   link->next       = DMForestTypeList;
6827d4645fSToby Isaac   DMForestTypeList = link;
6927d4645fSToby Isaac   PetscFunctionReturn(0);
7027d4645fSToby Isaac }
7127d4645fSToby Isaac 
729be51f97SToby Isaac /*@
739be51f97SToby Isaac   DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes
749be51f97SToby Isaac 
759be51f97SToby Isaac   Not Collective
769be51f97SToby Isaac 
779be51f97SToby Isaac   Input parameter:
789be51f97SToby Isaac . dm - the DM object
799be51f97SToby Isaac 
809be51f97SToby Isaac   Output parameter:
819be51f97SToby Isaac . isForest - whether dm is a subtype of DMFOREST
829be51f97SToby Isaac 
839be51f97SToby Isaac   Level: intermediate
849be51f97SToby Isaac 
859be51f97SToby Isaac .seealso: DMFOREST, DMForestRegisterType()
869be51f97SToby Isaac @*/
8727d4645fSToby Isaac PetscErrorCode DMIsForest(DM dm, PetscBool *isForest)
8827d4645fSToby Isaac {
8927d4645fSToby Isaac   DMForestTypeLink link = DMForestTypeList;
9027d4645fSToby Isaac   PetscErrorCode   ierr;
9127d4645fSToby Isaac 
9227d4645fSToby Isaac   PetscFunctionBegin;
9327d4645fSToby Isaac   while (link) {
9427d4645fSToby Isaac     PetscBool sameType;
9527d4645fSToby Isaac     ierr = PetscObjectTypeCompare((PetscObject)dm,link->name,&sameType);CHKERRQ(ierr);
9627d4645fSToby Isaac     if (sameType) {
9727d4645fSToby Isaac       *isForest = PETSC_TRUE;
9827d4645fSToby Isaac       PetscFunctionReturn(0);
9927d4645fSToby Isaac     }
10027d4645fSToby Isaac     link = link->next;
10127d4645fSToby Isaac   }
10227d4645fSToby Isaac   *isForest = PETSC_FALSE;
10327d4645fSToby Isaac   PetscFunctionReturn(0);
10427d4645fSToby Isaac }
10527d4645fSToby Isaac 
1069be51f97SToby Isaac /*@
1079be51f97SToby Isaac   DMForestTemplate - Create a new DM that will be adapted from a source DM.  The new DM reproduces the configuration
1089be51f97SToby 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
1099be51f97SToby Isaac   (by, e.g., refinement or repartitioning).  The source DM is also set as the adaptivity source DM of the new DM (see
1109be51f97SToby Isaac   DMForestSetAdaptivityForest()).
1119be51f97SToby Isaac 
1129be51f97SToby Isaac   Collective on dm
1139be51f97SToby Isaac 
1149be51f97SToby Isaac   Input Parameters:
1159be51f97SToby Isaac + dm - the source DM object
1169be51f97SToby Isaac - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen())
1179be51f97SToby Isaac 
1189be51f97SToby Isaac   Output Parameter:
1199be51f97SToby Isaac . tdm - the new DM object
1209be51f97SToby Isaac 
1219be51f97SToby Isaac   Level: intermediate
1229be51f97SToby Isaac 
1239be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest()
1249be51f97SToby Isaac @*/
1259be51f97SToby Isaac PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm)
126a0452a8eSToby Isaac {
127a0452a8eSToby Isaac   DM_Forest                  *forest = (DM_Forest*) dm->data;
12820e8089bSToby Isaac   DMType                     type;
129a0452a8eSToby Isaac   DM                         base;
130a0452a8eSToby Isaac   DMForestTopology           topology;
131a0452a8eSToby Isaac   PetscInt                   dim, overlap, ref, factor;
132a0452a8eSToby Isaac   DMForestAdaptivityStrategy strat;
133795844e7SToby Isaac   PetscDS                    ds;
134795844e7SToby Isaac   void                       *ctx;
13549fc9a2fSToby Isaac   PetscErrorCode             (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void*);
1363e58adeeSToby Isaac   void                       *mapCtx;
137a0452a8eSToby Isaac   PetscErrorCode             ierr;
138a0452a8eSToby Isaac 
139a0452a8eSToby Isaac   PetscFunctionBegin;
140a0452a8eSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14120e8089bSToby Isaac   ierr = DMCreate(PetscObjectComm((PetscObject)dm),tdm);CHKERRQ(ierr);
14220e8089bSToby Isaac   ierr = DMGetType(dm,&type);CHKERRQ(ierr);
14320e8089bSToby Isaac   ierr = DMSetType(*tdm,type);CHKERRQ(ierr);
144a0452a8eSToby Isaac   ierr = DMForestGetBaseDM(dm,&base);CHKERRQ(ierr);
14520e8089bSToby Isaac   ierr = DMForestSetBaseDM(*tdm,base);CHKERRQ(ierr);
146a0452a8eSToby Isaac   ierr = DMForestGetTopology(dm,&topology);CHKERRQ(ierr);
14720e8089bSToby Isaac   ierr = DMForestSetTopology(*tdm,topology);CHKERRQ(ierr);
148a0452a8eSToby Isaac   ierr = DMForestGetAdjacencyDimension(dm,&dim);CHKERRQ(ierr);
14920e8089bSToby Isaac   ierr = DMForestSetAdjacencyDimension(*tdm,dim);CHKERRQ(ierr);
150a0452a8eSToby Isaac   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
15120e8089bSToby Isaac   ierr = DMForestSetPartitionOverlap(*tdm,overlap);CHKERRQ(ierr);
152a0452a8eSToby Isaac   ierr = DMForestGetMinimumRefinement(dm,&ref);CHKERRQ(ierr);
15320e8089bSToby Isaac   ierr = DMForestSetMinimumRefinement(*tdm,ref);CHKERRQ(ierr);
154a0452a8eSToby Isaac   ierr = DMForestGetMaximumRefinement(dm,&ref);CHKERRQ(ierr);
15520e8089bSToby Isaac   ierr = DMForestSetMaximumRefinement(*tdm,ref);CHKERRQ(ierr);
156a0452a8eSToby Isaac   ierr = DMForestGetAdaptivityStrategy(dm,&strat);CHKERRQ(ierr);
15720e8089bSToby Isaac   ierr = DMForestSetAdaptivityStrategy(*tdm,strat);CHKERRQ(ierr);
158a0452a8eSToby Isaac   ierr = DMForestGetGradeFactor(dm,&factor);CHKERRQ(ierr);
15920e8089bSToby Isaac   ierr = DMForestSetGradeFactor(*tdm,factor);CHKERRQ(ierr);
1603e58adeeSToby Isaac   ierr = DMForestGetBaseCoordinateMapping(dm,&map,&mapCtx);CHKERRQ(ierr);
1615e8c540aSToby Isaac   ierr = DMForestSetBaseCoordinateMapping(*tdm,map,mapCtx);CHKERRQ(ierr);
162a0452a8eSToby Isaac   if (forest->ftemplate) {
16320e8089bSToby Isaac     ierr = (forest->ftemplate)(dm, *tdm);CHKERRQ(ierr);
164a0452a8eSToby Isaac   }
16520e8089bSToby Isaac   ierr = DMForestSetAdaptivityForest(*tdm,dm);CHKERRQ(ierr);
166795844e7SToby Isaac   ierr = DMGetDS(dm,&ds);CHKERRQ(ierr);
167795844e7SToby Isaac   ierr = DMSetDS(*tdm,ds);CHKERRQ(ierr);
168795844e7SToby Isaac   ierr = DMGetApplicationContext(dm,&ctx);CHKERRQ(ierr);
169795844e7SToby Isaac   ierr = DMSetApplicationContext(*tdm,&ctx);CHKERRQ(ierr);
17090b157c4SStefano Zampini   {
17190b157c4SStefano Zampini     PetscBool            isper;
172795844e7SToby Isaac     const PetscReal      *maxCell, *L;
173795844e7SToby Isaac     const DMBoundaryType *bd;
174795844e7SToby Isaac 
17590b157c4SStefano Zampini     ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
17690b157c4SStefano Zampini     ierr = DMSetPeriodicity(*tdm,isper,maxCell,L,bd);CHKERRQ(ierr);
177795844e7SToby Isaac   }
178bff67a9bSToby Isaac   ierr = DMCopyBoundary(dm,*tdm);CHKERRQ(ierr);
179a0452a8eSToby Isaac   PetscFunctionReturn(0);
180a0452a8eSToby Isaac }
181a0452a8eSToby Isaac 
18201d9d024SToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm);
18301d9d024SToby Isaac 
184db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm)
185db4d5e8cSToby Isaac {
186db4d5e8cSToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
187db4d5e8cSToby Isaac   const char     *type;
188db4d5e8cSToby Isaac   PetscErrorCode ierr;
189db4d5e8cSToby Isaac 
190db4d5e8cSToby Isaac   PetscFunctionBegin;
191db4d5e8cSToby Isaac   forest->refct++;
192db4d5e8cSToby Isaac   (*newdm)->data = forest;
193db4d5e8cSToby Isaac   ierr           = PetscObjectGetType((PetscObject) dm, &type);CHKERRQ(ierr);
194db4d5e8cSToby Isaac   ierr           = PetscObjectChangeTypeName((PetscObject) *newdm, type);CHKERRQ(ierr);
19501d9d024SToby Isaac   ierr           = DMInitialize_Forest(*newdm);CHKERRQ(ierr);
196db4d5e8cSToby Isaac   PetscFunctionReturn(0);
197db4d5e8cSToby Isaac }
198db4d5e8cSToby Isaac 
199d222f98bSToby Isaac static PetscErrorCode DMDestroy_Forest(DM dm)
200db4d5e8cSToby Isaac {
201db4d5e8cSToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
202db4d5e8cSToby Isaac   PetscErrorCode ierr;
203db4d5e8cSToby Isaac 
204db4d5e8cSToby Isaac   PetscFunctionBegin;
205db4d5e8cSToby Isaac   if (--forest->refct > 0) PetscFunctionReturn(0);
206d222f98bSToby Isaac   if (forest->destroy) {ierr = forest->destroy(dm);CHKERRQ(ierr);}
207db4d5e8cSToby Isaac   ierr = PetscSFDestroy(&forest->cellSF);CHKERRQ(ierr);
2080f17b9e3SToby Isaac   ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr);
2090f17b9e3SToby Isaac   ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr);
210a1b0c543SToby Isaac   ierr = DMLabelDestroy(&forest->adaptLabel);CHKERRQ(ierr);
2119a81d013SToby Isaac   ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr);
21256ba9f64SToby Isaac   ierr = DMDestroy(&forest->base);CHKERRQ(ierr);
213ba936b91SToby Isaac   ierr = DMDestroy(&forest->adapt);CHKERRQ(ierr);
21430f902e7SToby Isaac   ierr = PetscFree(forest->topology);CHKERRQ(ierr);
21530f902e7SToby Isaac   ierr = PetscFree(forest);CHKERRQ(ierr);
216db4d5e8cSToby Isaac   PetscFunctionReturn(0);
217db4d5e8cSToby Isaac }
218db4d5e8cSToby Isaac 
2199be51f97SToby Isaac /*@C
2209be51f97SToby Isaac   DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase.  The topology is a string (e.g.
22168d54884SBarry Smith   "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest during
2229be51f97SToby Isaac   DMSetUp().
2239be51f97SToby Isaac 
2249be51f97SToby Isaac   Logically collective on dm
2259be51f97SToby Isaac 
2269be51f97SToby Isaac   Input parameters:
2279be51f97SToby Isaac + dm - the forest
2289be51f97SToby Isaac - topology - the topology of the forest
2299be51f97SToby Isaac 
2309be51f97SToby Isaac   Level: intermediate
2319be51f97SToby Isaac 
2329be51f97SToby Isaac .seealso(): DMForestGetTopology(), DMForestSetBaseDM()
2339be51f97SToby Isaac @*/
234dd8e54a2SToby Isaac PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology)
235db4d5e8cSToby Isaac {
236db4d5e8cSToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
237db4d5e8cSToby Isaac   PetscErrorCode ierr;
238db4d5e8cSToby Isaac 
239db4d5e8cSToby Isaac   PetscFunctionBegin;
240db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
241ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup");
242dd8e54a2SToby Isaac   ierr = PetscFree(forest->topology);CHKERRQ(ierr);
243dd8e54a2SToby Isaac   ierr = PetscStrallocpy((const char*)topology,(char**) &forest->topology);CHKERRQ(ierr);
244db4d5e8cSToby Isaac   PetscFunctionReturn(0);
245db4d5e8cSToby Isaac }
246db4d5e8cSToby Isaac 
2479be51f97SToby Isaac /*@C
2489be51f97SToby Isaac   DMForestGetTopology - Get a string describing the topology of a DMForest.
2499be51f97SToby Isaac 
2509be51f97SToby Isaac   Not collective
2519be51f97SToby Isaac 
2529be51f97SToby Isaac   Input parameter:
2539be51f97SToby Isaac . dm - the forest
2549be51f97SToby Isaac 
2559be51f97SToby Isaac   Output parameter:
2569be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell')
2579be51f97SToby Isaac 
2589be51f97SToby Isaac   Level: intermediate
2599be51f97SToby Isaac 
2609be51f97SToby Isaac .seealso: DMForestSetTopology()
2619be51f97SToby Isaac @*/
262dd8e54a2SToby Isaac PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology)
263dd8e54a2SToby Isaac {
264dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
265dd8e54a2SToby Isaac 
266dd8e54a2SToby Isaac   PetscFunctionBegin;
267dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
268dd8e54a2SToby Isaac   PetscValidPointer(topology,2);
269dd8e54a2SToby Isaac   *topology = forest->topology;
270dd8e54a2SToby Isaac   PetscFunctionReturn(0);
271dd8e54a2SToby Isaac }
272dd8e54a2SToby Isaac 
2739be51f97SToby Isaac /*@
2749be51f97SToby Isaac   DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest.  The
2759be51f97SToby Isaac   forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its
276765b024eSBarry Smith   base.  In general, two forest must share a base to be comparable, to do things like construct interpolators.
2779be51f97SToby Isaac 
2789be51f97SToby Isaac   Logically collective on dm
2799be51f97SToby Isaac 
2809be51f97SToby Isaac   Input Parameters:
2819be51f97SToby Isaac + dm - the forest
2829be51f97SToby Isaac - base - the base DM of the forest
2839be51f97SToby Isaac 
284765b024eSBarry Smith   Notes: Currently the base DM must be a DMPLEX
285765b024eSBarry Smith 
2869be51f97SToby Isaac   Level: intermediate
2879be51f97SToby Isaac 
2889be51f97SToby Isaac .seealso(): DMForestGetBaseDM()
2899be51f97SToby Isaac @*/
290dd8e54a2SToby Isaac PetscErrorCode DMForestSetBaseDM(DM dm, DM base)
291dd8e54a2SToby Isaac {
292dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
293dd8e54a2SToby Isaac   PetscInt       dim, dimEmbed;
294dd8e54a2SToby Isaac   PetscErrorCode ierr;
295dd8e54a2SToby Isaac 
296dd8e54a2SToby Isaac   PetscFunctionBegin;
297dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
298ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup");
299dd8e54a2SToby Isaac   ierr         = PetscObjectReference((PetscObject)base);CHKERRQ(ierr);
300dd8e54a2SToby Isaac   ierr         = DMDestroy(&forest->base);CHKERRQ(ierr);
301dd8e54a2SToby Isaac   forest->base = base;
302a0452a8eSToby Isaac   if (base) {
303a0452a8eSToby Isaac     PetscValidHeaderSpecific(base, DM_CLASSID, 2);
304dd8e54a2SToby Isaac     ierr = DMGetDimension(base,&dim);CHKERRQ(ierr);
305dd8e54a2SToby Isaac     ierr = DMSetDimension(dm,dim);CHKERRQ(ierr);
306dd8e54a2SToby Isaac     ierr = DMGetCoordinateDim(base,&dimEmbed);CHKERRQ(ierr);
307dd8e54a2SToby Isaac     ierr = DMSetCoordinateDim(dm,dimEmbed);CHKERRQ(ierr);
308a0452a8eSToby Isaac   }
309dd8e54a2SToby Isaac   PetscFunctionReturn(0);
310dd8e54a2SToby Isaac }
311dd8e54a2SToby Isaac 
3129be51f97SToby Isaac /*@
3139be51f97SToby Isaac   DMForestGetBaseDM - Get the base DM of a DMForest forest.  The forest will be hierarchically refined from the base,
31468d54884SBarry Smith   and all refinements/coarsenings of the forest will share its base.  In general, two forest must share a base to be
3159be51f97SToby Isaac   comparable, to do things like construct interpolators.
3169be51f97SToby Isaac 
3179be51f97SToby Isaac   Not collective
3189be51f97SToby Isaac 
3199be51f97SToby Isaac   Input Parameter:
3209be51f97SToby Isaac . dm - the forest
3219be51f97SToby Isaac 
3229be51f97SToby Isaac   Output Parameter:
3239be51f97SToby Isaac . base - the base DM of the forest
3249be51f97SToby Isaac 
3259be51f97SToby Isaac   Level: intermediate
3269be51f97SToby Isaac 
3279be51f97SToby Isaac .seealso(); DMForestSetBaseDM()
3289be51f97SToby Isaac @*/
329dd8e54a2SToby Isaac PetscErrorCode DMForestGetBaseDM(DM dm, DM *base)
330dd8e54a2SToby Isaac {
331dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
332dd8e54a2SToby Isaac 
333dd8e54a2SToby Isaac   PetscFunctionBegin;
334dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
335dd8e54a2SToby Isaac   PetscValidPointer(base, 2);
336dd8e54a2SToby Isaac   *base = forest->base;
337dd8e54a2SToby Isaac   PetscFunctionReturn(0);
338dd8e54a2SToby Isaac }
339dd8e54a2SToby Isaac 
34099478f86SToby Isaac PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
341cf38a08cSToby Isaac {
342cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
343cf38a08cSToby Isaac 
344cf38a08cSToby Isaac   PetscFunctionBegin;
345cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
346cf38a08cSToby Isaac   forest->mapcoordinates    = func;
347cf38a08cSToby Isaac   forest->mapcoordinatesctx = ctx;
348cf38a08cSToby Isaac   PetscFunctionReturn(0);
349cf38a08cSToby Isaac }
350cf38a08cSToby Isaac 
35199478f86SToby Isaac PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func) (DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
352cf38a08cSToby Isaac {
353cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
354cf38a08cSToby Isaac 
355cf38a08cSToby Isaac   PetscFunctionBegin;
356cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
357cf38a08cSToby Isaac   if (func) *func = forest->mapcoordinates;
358cf38a08cSToby Isaac   if (ctx) *((void**) ctx) = forest->mapcoordinatesctx;
359cf38a08cSToby Isaac   PetscFunctionReturn(0);
360cf38a08cSToby Isaac }
361cf38a08cSToby Isaac 
3629be51f97SToby Isaac /*@
3639be51f97SToby Isaac   DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be
3649be51f97SToby Isaac   adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp().  Usually not needed
3659be51f97SToby Isaac   by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this
3669be51f97SToby Isaac   routine.
3679be51f97SToby Isaac 
368dffe73a3SToby Isaac   Note that this can be called after setup with adapt = NULL, which will clear all internal data related to the
369dffe73a3SToby Isaac   adaptivity forest from dm.  This way, repeatedly adapting does not leave stale DM objects in memory.
370dffe73a3SToby Isaac 
3719be51f97SToby Isaac   Logically collective on dm
3729be51f97SToby Isaac 
3739be51f97SToby Isaac   Input Parameter:
3749be51f97SToby Isaac + dm - the new forest, which will be constructed from adapt
3759be51f97SToby Isaac - adapt - the old forest
3769be51f97SToby Isaac 
3779be51f97SToby Isaac   Level: intermediate
3789be51f97SToby Isaac 
3799be51f97SToby Isaac .seealso: DMForestGetAdaptivityForest(), DMForestSetAdaptivityPurpose()
3809be51f97SToby Isaac @*/
381ba936b91SToby Isaac PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt)
382dd8e54a2SToby Isaac {
383dffe73a3SToby Isaac   DM_Forest      *forest, *adaptForest, *oldAdaptForest;
384dffe73a3SToby Isaac   DM             oldAdapt;
385*456cc5b7SMatthew G. Knepley   PetscBool      isForest;
386dd8e54a2SToby Isaac   PetscErrorCode ierr;
387dd8e54a2SToby Isaac 
388dd8e54a2SToby Isaac   PetscFunctionBegin;
389dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
390ba936b91SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
391*456cc5b7SMatthew G. Knepley   ierr = DMIsForest(dm, &isForest);CHKERRQ(ierr);
392*456cc5b7SMatthew G. Knepley   if (!isForest) PetscFunctionReturn(0);
393ba936b91SToby Isaac   forest   = (DM_Forest*) dm->data;
394dffe73a3SToby Isaac   ierr     = DMForestGetAdaptivityForest(dm,&oldAdapt);CHKERRQ(ierr);
395dffe73a3SToby Isaac   if (adapt != NULL && dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
396193eb951SToby Isaac   adaptForest    = (DM_Forest*) (adapt ? adapt->data : NULL);
397193eb951SToby Isaac   oldAdaptForest = (DM_Forest*) (oldAdapt ? oldAdapt->data : NULL);
398dffe73a3SToby Isaac   if (adaptForest != oldAdaptForest) {
399dffe73a3SToby Isaac     ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr);
400dffe73a3SToby Isaac     ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr);
401dffe73a3SToby Isaac     if (forest->clearadaptivityforest) {ierr = (forest->clearadaptivityforest)(dm);CHKERRQ(ierr);}
402dffe73a3SToby Isaac   }
40326d9498aSToby Isaac   switch (forest->adaptPurpose) {
404cd3c525cSToby Isaac   case DM_ADAPT_DETERMINE:
405ba936b91SToby Isaac     ierr          = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
406ba936b91SToby Isaac     ierr          = DMDestroy(&(forest->adapt));CHKERRQ(ierr);
407ba936b91SToby Isaac     forest->adapt = adapt;
40826d9498aSToby Isaac     break;
409a1b0c543SToby Isaac   case DM_ADAPT_REFINE:
41026d9498aSToby Isaac     ierr = DMSetCoarseDM(dm,adapt);CHKERRQ(ierr);
41126d9498aSToby Isaac     break;
412a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
41326d9498aSToby Isaac     ierr = DMSetFineDM(dm,adapt);CHKERRQ(ierr);
41426d9498aSToby Isaac     break;
41526d9498aSToby Isaac   default:
41626d9498aSToby Isaac     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
41726d9498aSToby Isaac   }
418dd8e54a2SToby Isaac   PetscFunctionReturn(0);
419dd8e54a2SToby Isaac }
420dd8e54a2SToby Isaac 
4219be51f97SToby Isaac /*@
4229be51f97SToby Isaac   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.
4239be51f97SToby Isaac 
4249be51f97SToby Isaac   Not collective
4259be51f97SToby Isaac 
4269be51f97SToby Isaac   Input Parameter:
4279be51f97SToby Isaac . dm - the forest
4289be51f97SToby Isaac 
4299be51f97SToby Isaac   Output Parameter:
4309be51f97SToby Isaac . adapt - the forest from which dm is/was adapted
4319be51f97SToby Isaac 
4329be51f97SToby Isaac   Level: intermediate
4339be51f97SToby Isaac 
4349be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose()
4359be51f97SToby Isaac @*/
436ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
437dd8e54a2SToby Isaac {
438ba936b91SToby Isaac   DM_Forest      *forest;
43926d9498aSToby Isaac   PetscErrorCode ierr;
440dd8e54a2SToby Isaac 
441dd8e54a2SToby Isaac   PetscFunctionBegin;
442dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
443ba936b91SToby Isaac   forest = (DM_Forest*) dm->data;
44426d9498aSToby Isaac   switch (forest->adaptPurpose) {
445cd3c525cSToby Isaac   case DM_ADAPT_DETERMINE:
446ba936b91SToby Isaac     *adapt = forest->adapt;
44726d9498aSToby Isaac     break;
448a1b0c543SToby Isaac   case DM_ADAPT_REFINE:
44926d9498aSToby Isaac     ierr = DMGetCoarseDM(dm,adapt);CHKERRQ(ierr);
45026d9498aSToby Isaac     break;
451a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
45226d9498aSToby Isaac     ierr = DMGetFineDM(dm,adapt);CHKERRQ(ierr);
45326d9498aSToby Isaac     break;
45426d9498aSToby Isaac   default:
45526d9498aSToby Isaac     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
45626d9498aSToby Isaac   }
45726d9498aSToby Isaac   PetscFunctionReturn(0);
45826d9498aSToby Isaac }
45926d9498aSToby Isaac 
4609be51f97SToby Isaac /*@
4619be51f97SToby Isaac   DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its
462a1b0c543SToby Isaac   source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening
463cd3c525cSToby Isaac   (DM_ADAPT_COARSEN), or undefined (DM_ADAPT_DETERMINE).  This only matters for the purposes of reference counting:
4649be51f97SToby Isaac   during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse
4659be51f97SToby Isaac   relationship (see DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does
4669be51f97SToby Isaac   not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can
4679be51f97SToby Isaac   cause a memory leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
4689be51f97SToby Isaac 
4699be51f97SToby Isaac   Logically collective on dm
4709be51f97SToby Isaac 
4719be51f97SToby Isaac   Input Parameters:
4729be51f97SToby Isaac + dm - the forest
473cd3c525cSToby Isaac - purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN)
4749be51f97SToby Isaac 
4759be51f97SToby Isaac   Level: advanced
4769be51f97SToby Isaac 
4779be51f97SToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
4789be51f97SToby Isaac @*/
479a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose)
48026d9498aSToby Isaac {
48126d9498aSToby Isaac   DM_Forest      *forest;
48226d9498aSToby Isaac   PetscErrorCode ierr;
48326d9498aSToby Isaac 
48426d9498aSToby Isaac   PetscFunctionBegin;
48526d9498aSToby Isaac   forest = (DM_Forest*) dm->data;
4869be51f97SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
48726d9498aSToby Isaac   if (purpose != forest->adaptPurpose) {
48826d9498aSToby Isaac     DM adapt;
48926d9498aSToby Isaac 
49026d9498aSToby Isaac     ierr = DMForestGetAdaptivityForest(dm,&adapt);CHKERRQ(ierr);
49126d9498aSToby Isaac     ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
49226d9498aSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
493f885a11aSToby Isaac 
49426d9498aSToby Isaac     forest->adaptPurpose = purpose;
495f885a11aSToby Isaac 
49626d9498aSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,adapt);CHKERRQ(ierr);
49726d9498aSToby Isaac     ierr = DMDestroy(&adapt);CHKERRQ(ierr);
49826d9498aSToby Isaac   }
499dd8e54a2SToby Isaac   PetscFunctionReturn(0);
500dd8e54a2SToby Isaac }
501dd8e54a2SToby Isaac 
50256c0450aSToby Isaac /*@
50356c0450aSToby Isaac   DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with
504a1b0c543SToby Isaac   DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening (DM_ADAPT_COARSEN), or
505cd3c525cSToby Isaac   undefined (DM_ADAPT_DETERMINE).  This only matters for the purposes of reference counting: during DMDestroy(), cyclic
50656c0450aSToby Isaac   references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see
50756c0450aSToby Isaac   DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does not maintain a
50856c0450aSToby Isaac   reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory
50956c0450aSToby Isaac   leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
51056c0450aSToby Isaac 
51156c0450aSToby Isaac   Not collective
51256c0450aSToby Isaac 
51356c0450aSToby Isaac   Input Parameter:
51456c0450aSToby Isaac . dm - the forest
51556c0450aSToby Isaac 
51656c0450aSToby Isaac   Output Parameter:
517cd3c525cSToby Isaac . purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN)
51856c0450aSToby Isaac 
51956c0450aSToby Isaac   Level: advanced
52056c0450aSToby Isaac 
52156c0450aSToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
52256c0450aSToby Isaac @*/
523a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose)
52456c0450aSToby Isaac {
52556c0450aSToby Isaac   DM_Forest *forest;
52656c0450aSToby Isaac 
52756c0450aSToby Isaac   PetscFunctionBegin;
52856c0450aSToby Isaac   forest   = (DM_Forest*) dm->data;
52956c0450aSToby Isaac   *purpose = forest->adaptPurpose;
53056c0450aSToby Isaac   PetscFunctionReturn(0);
53156c0450aSToby Isaac }
53256c0450aSToby Isaac 
5339be51f97SToby Isaac /*@
5349be51f97SToby Isaac   DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
5359be51f97SToby Isaac   cell adjacency (for the purposes of partitioning and overlap).
5369be51f97SToby Isaac 
5379be51f97SToby Isaac   Logically collective on dm
5389be51f97SToby Isaac 
5399be51f97SToby Isaac   Input Parameters:
5409be51f97SToby Isaac + dm - the forest
5419be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency)
5429be51f97SToby Isaac 
5439be51f97SToby Isaac   Level: intermediate
5449be51f97SToby Isaac 
5459be51f97SToby Isaac .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap()
5469be51f97SToby Isaac @*/
547dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
548dd8e54a2SToby Isaac {
549dd8e54a2SToby Isaac   PetscInt       dim;
550dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
551dd8e54a2SToby Isaac   PetscErrorCode ierr;
552dd8e54a2SToby Isaac 
553dd8e54a2SToby Isaac   PetscFunctionBegin;
554dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
555ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup");
556dd8e54a2SToby Isaac   if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim);
557dd8e54a2SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
558dd8e54a2SToby Isaac   if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim);
559dd8e54a2SToby Isaac   forest->adjDim = adjDim;
560dd8e54a2SToby Isaac   PetscFunctionReturn(0);
561dd8e54a2SToby Isaac }
562dd8e54a2SToby Isaac 
5639be51f97SToby Isaac /*@
5649be51f97SToby Isaac   DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that,
5659be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
5669be51f97SToby Isaac 
5679be51f97SToby Isaac   Logically collective on dm
5689be51f97SToby Isaac 
5699be51f97SToby Isaac   Input Parameters:
5709be51f97SToby Isaac + dm - the forest
5719be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
5729be51f97SToby Isaac 
5739be51f97SToby Isaac   Level: intermediate
5749be51f97SToby Isaac 
5759be51f97SToby Isaac .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension()
5769be51f97SToby Isaac @*/
577dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
578dd8e54a2SToby Isaac {
579dd8e54a2SToby Isaac   PetscInt       dim;
580dd8e54a2SToby Isaac   PetscErrorCode ierr;
581dd8e54a2SToby Isaac 
582dd8e54a2SToby Isaac   PetscFunctionBegin;
583dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
584dd8e54a2SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
585dd8e54a2SToby Isaac   ierr = DMForestSetAdjacencyDimension(dm,dim-adjCodim);CHKERRQ(ierr);
586dd8e54a2SToby Isaac   PetscFunctionReturn(0);
587dd8e54a2SToby Isaac }
588dd8e54a2SToby Isaac 
5899be51f97SToby Isaac /*@
5909be51f97SToby Isaac   DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
5919be51f97SToby Isaac   purposes of partitioning and overlap).
5929be51f97SToby Isaac 
5939be51f97SToby Isaac   Not collective
5949be51f97SToby Isaac 
5959be51f97SToby Isaac   Input Parameter:
5969be51f97SToby Isaac . dm - the forest
5979be51f97SToby Isaac 
5989be51f97SToby Isaac   Output Parameter:
5999be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency)
6009be51f97SToby Isaac 
6019be51f97SToby Isaac   Level: intermediate
6029be51f97SToby Isaac 
6039be51f97SToby Isaac .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap()
6049be51f97SToby Isaac @*/
605dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
606dd8e54a2SToby Isaac {
607dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
608dd8e54a2SToby Isaac 
609dd8e54a2SToby Isaac   PetscFunctionBegin;
610dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
611dd8e54a2SToby Isaac   PetscValidIntPointer(adjDim,2);
612dd8e54a2SToby Isaac   *adjDim = forest->adjDim;
613dd8e54a2SToby Isaac   PetscFunctionReturn(0);
614dd8e54a2SToby Isaac }
615dd8e54a2SToby Isaac 
6169be51f97SToby Isaac /*@
6179be51f97SToby Isaac   DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that,
6189be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
6199be51f97SToby Isaac 
6209be51f97SToby Isaac   Not collective
6219be51f97SToby Isaac 
6229be51f97SToby Isaac   Input Parameter:
6239be51f97SToby Isaac . dm - the forest
6249be51f97SToby Isaac 
6259be51f97SToby Isaac   Output Parameter:
6269be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
6279be51f97SToby Isaac 
6289be51f97SToby Isaac   Level: intermediate
6299be51f97SToby Isaac 
6309be51f97SToby Isaac .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension()
6319be51f97SToby Isaac @*/
632dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
633dd8e54a2SToby Isaac {
634dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
635dd8e54a2SToby Isaac   PetscInt       dim;
636dd8e54a2SToby Isaac   PetscErrorCode ierr;
637dd8e54a2SToby Isaac 
638dd8e54a2SToby Isaac   PetscFunctionBegin;
639dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
640dd8e54a2SToby Isaac   PetscValidIntPointer(adjCodim,2);
641dd8e54a2SToby Isaac   ierr      = DMGetDimension(dm,&dim);CHKERRQ(ierr);
642dd8e54a2SToby Isaac   *adjCodim = dim - forest->adjDim;
643dd8e54a2SToby Isaac   PetscFunctionReturn(0);
644dd8e54a2SToby Isaac }
645dd8e54a2SToby Isaac 
6469be51f97SToby Isaac /*@
6479be51f97SToby Isaac   DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel
6489be51f97SToby Isaac   partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding
6499be51f97SToby Isaac   adjacent cells
6509be51f97SToby Isaac 
6519be51f97SToby Isaac   Logically collective on dm
6529be51f97SToby Isaac 
6539be51f97SToby Isaac   Input Parameters:
6549be51f97SToby Isaac + dm - the forest
6559be51f97SToby Isaac - overlap - default 0
6569be51f97SToby Isaac 
6579be51f97SToby Isaac   Level: intermediate
6589be51f97SToby Isaac 
6599be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
6609be51f97SToby Isaac @*/
661dd8e54a2SToby Isaac PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap)
662dd8e54a2SToby Isaac {
663dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
664dd8e54a2SToby Isaac 
665dd8e54a2SToby Isaac   PetscFunctionBegin;
666dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
667ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup");
668dd8e54a2SToby Isaac   if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap);
669dd8e54a2SToby Isaac   forest->overlap = overlap;
670dd8e54a2SToby Isaac   PetscFunctionReturn(0);
671dd8e54a2SToby Isaac }
672dd8e54a2SToby Isaac 
6739be51f97SToby Isaac /*@
6749be51f97SToby Isaac   DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values
6759be51f97SToby Isaac   > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells
6769be51f97SToby Isaac 
6779be51f97SToby Isaac   Not collective
6789be51f97SToby Isaac 
6799be51f97SToby Isaac   Input Parameter:
6809be51f97SToby Isaac . dm - the forest
6819be51f97SToby Isaac 
6829be51f97SToby Isaac   Output Parameter:
6839be51f97SToby Isaac . overlap - default 0
6849be51f97SToby Isaac 
6859be51f97SToby Isaac   Level: intermediate
6869be51f97SToby Isaac 
6879be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
6889be51f97SToby Isaac @*/
689dd8e54a2SToby Isaac PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap)
690dd8e54a2SToby Isaac {
691dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
692dd8e54a2SToby Isaac 
693dd8e54a2SToby Isaac   PetscFunctionBegin;
694dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
695dd8e54a2SToby Isaac   PetscValidIntPointer(overlap,2);
696dd8e54a2SToby Isaac   *overlap = forest->overlap;
697dd8e54a2SToby Isaac   PetscFunctionReturn(0);
698dd8e54a2SToby Isaac }
699dd8e54a2SToby Isaac 
7009be51f97SToby Isaac /*@
7019be51f97SToby Isaac   DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base
7029be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest
7039be51f97SToby Isaac   (see DMForestGetAdaptivityForest()) this limits the amount of coarsening.
7049be51f97SToby Isaac 
7059be51f97SToby Isaac   Logically collective on dm
7069be51f97SToby Isaac 
7079be51f97SToby Isaac   Input Parameters:
7089be51f97SToby Isaac + dm - the forest
7099be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7109be51f97SToby Isaac 
7119be51f97SToby Isaac   Level: intermediate
7129be51f97SToby Isaac 
7139be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
7149be51f97SToby Isaac @*/
715dd8e54a2SToby Isaac PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement)
716dd8e54a2SToby Isaac {
717dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
718dd8e54a2SToby Isaac 
719dd8e54a2SToby Isaac   PetscFunctionBegin;
720dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
721ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup");
722dd8e54a2SToby Isaac   forest->minRefinement = minRefinement;
723dd8e54a2SToby Isaac   PetscFunctionReturn(0);
724dd8e54a2SToby Isaac }
725dd8e54a2SToby Isaac 
7269be51f97SToby Isaac /*@
7279be51f97SToby Isaac   DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see
7289be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest (see
7299be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of coarsening.
7309be51f97SToby Isaac 
7319be51f97SToby Isaac   Not collective
7329be51f97SToby Isaac 
7339be51f97SToby Isaac   Input Parameter:
7349be51f97SToby Isaac . dm - the forest
7359be51f97SToby Isaac 
7369be51f97SToby Isaac   Output Parameter:
7379be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7389be51f97SToby Isaac 
7399be51f97SToby Isaac   Level: intermediate
7409be51f97SToby Isaac 
7419be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
7429be51f97SToby Isaac @*/
743dd8e54a2SToby Isaac PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement)
744dd8e54a2SToby Isaac {
745dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
746dd8e54a2SToby Isaac 
747dd8e54a2SToby Isaac   PetscFunctionBegin;
748dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
749dd8e54a2SToby Isaac   PetscValidIntPointer(minRefinement,2);
750dd8e54a2SToby Isaac   *minRefinement = forest->minRefinement;
751dd8e54a2SToby Isaac   PetscFunctionReturn(0);
752dd8e54a2SToby Isaac }
753dd8e54a2SToby Isaac 
7549be51f97SToby Isaac /*@
7559be51f97SToby Isaac   DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base
7569be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.
7579be51f97SToby Isaac 
7589be51f97SToby Isaac   Logically collective on dm
7599be51f97SToby Isaac 
7609be51f97SToby Isaac   Input Parameters:
7619be51f97SToby Isaac + dm - the forest
7629be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7639be51f97SToby Isaac 
7649be51f97SToby Isaac   Level: intermediate
7659be51f97SToby Isaac 
7669be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
7679be51f97SToby Isaac @*/
76856ba9f64SToby Isaac PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement)
76956ba9f64SToby Isaac {
77056ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
77156ba9f64SToby Isaac 
77256ba9f64SToby Isaac   PetscFunctionBegin;
77356ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
77456ba9f64SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup");
77556ba9f64SToby Isaac   forest->initRefinement = initRefinement;
77656ba9f64SToby Isaac   PetscFunctionReturn(0);
77756ba9f64SToby Isaac }
77856ba9f64SToby Isaac 
7799be51f97SToby Isaac /*@
7809be51f97SToby Isaac   DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see
7819be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.
7829be51f97SToby Isaac 
7839be51f97SToby Isaac   Not collective
7849be51f97SToby Isaac 
7859be51f97SToby Isaac   Input Parameter:
7869be51f97SToby Isaac . dm - the forest
7879be51f97SToby Isaac 
7889be51f97SToby Isaac   Output Paramater:
7899be51f97SToby Isaac . initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7909be51f97SToby Isaac 
7919be51f97SToby Isaac   Level: intermediate
7929be51f97SToby Isaac 
7939be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
7949be51f97SToby Isaac @*/
79556ba9f64SToby Isaac PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement)
79656ba9f64SToby Isaac {
79756ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
79856ba9f64SToby Isaac 
79956ba9f64SToby Isaac   PetscFunctionBegin;
80056ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
80156ba9f64SToby Isaac   PetscValidIntPointer(initRefinement,2);
80256ba9f64SToby Isaac   *initRefinement = forest->initRefinement;
80356ba9f64SToby Isaac   PetscFunctionReturn(0);
80456ba9f64SToby Isaac }
80556ba9f64SToby Isaac 
8069be51f97SToby Isaac /*@
8079be51f97SToby Isaac   DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base
8089be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest
8099be51f97SToby Isaac   (see DMForestGetAdaptivityForest()), this limits the amount of refinement.
8109be51f97SToby Isaac 
8119be51f97SToby Isaac   Logically collective on dm
8129be51f97SToby Isaac 
8139be51f97SToby Isaac   Input Parameters:
8149be51f97SToby Isaac + dm - the forest
8159be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8169be51f97SToby Isaac 
8179be51f97SToby Isaac   Level: intermediate
8189be51f97SToby Isaac 
8199be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM()
8209be51f97SToby Isaac @*/
821c7eeac06SToby Isaac PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement)
822dd8e54a2SToby Isaac {
823dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
824dd8e54a2SToby Isaac 
825dd8e54a2SToby Isaac   PetscFunctionBegin;
826dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
827ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup");
828c7eeac06SToby Isaac   forest->maxRefinement = maxRefinement;
829dd8e54a2SToby Isaac   PetscFunctionReturn(0);
830dd8e54a2SToby Isaac }
831dd8e54a2SToby Isaac 
8329be51f97SToby Isaac /*@
8339be51f97SToby Isaac   DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see
8349be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest (see
8359be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of refinement.
8369be51f97SToby Isaac 
8379be51f97SToby Isaac   Not collective
8389be51f97SToby Isaac 
8399be51f97SToby Isaac   Input Parameter:
8409be51f97SToby Isaac . dm - the forest
8419be51f97SToby Isaac 
8429be51f97SToby Isaac   Output Parameter:
8439be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8449be51f97SToby Isaac 
8459be51f97SToby Isaac   Level: intermediate
8469be51f97SToby Isaac 
8479be51f97SToby Isaac .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
8489be51f97SToby Isaac @*/
849c7eeac06SToby Isaac PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement)
850dd8e54a2SToby Isaac {
851dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
852dd8e54a2SToby Isaac 
853dd8e54a2SToby Isaac   PetscFunctionBegin;
854dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
855c7eeac06SToby Isaac   PetscValidIntPointer(maxRefinement,2);
856c7eeac06SToby Isaac   *maxRefinement = forest->maxRefinement;
857dd8e54a2SToby Isaac   PetscFunctionReturn(0);
858dd8e54a2SToby Isaac }
859c7eeac06SToby Isaac 
8609be51f97SToby Isaac /*@C
8619be51f97SToby Isaac   DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes.
8629be51f97SToby Isaac   Subtypes of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree
8639be51f97SToby Isaac   for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to
8649be51f97SToby Isaac   specify refinement/coarsening.
8659be51f97SToby Isaac 
8669be51f97SToby Isaac   Logically collective on dm
8679be51f97SToby Isaac 
8689be51f97SToby Isaac   Input Parameters:
8699be51f97SToby Isaac + dm - the forest
8709be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL
8719be51f97SToby Isaac 
8729be51f97SToby Isaac   Level: advanced
8739be51f97SToby Isaac 
8749be51f97SToby Isaac .seealso: DMForestGetAdaptivityStrategy()
8759be51f97SToby Isaac @*/
876c7eeac06SToby Isaac PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy)
877c7eeac06SToby Isaac {
878c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
879c7eeac06SToby Isaac   PetscErrorCode ierr;
880c7eeac06SToby Isaac 
881c7eeac06SToby Isaac   PetscFunctionBegin;
882c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
883c7eeac06SToby Isaac   ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr);
884a73e2921SToby Isaac   ierr = PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy);CHKERRQ(ierr);
885c7eeac06SToby Isaac   PetscFunctionReturn(0);
886c7eeac06SToby Isaac }
887c7eeac06SToby Isaac 
8889be51f97SToby Isaac /*@C
8899be51f97SToby Isaac   DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes.  Subtypes
8909be51f97SToby Isaac   of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all
8919be51f97SToby Isaac   processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only
8929be51f97SToby Isaac   one process needs to specify refinement/coarsening.
8939be51f97SToby Isaac 
8949be51f97SToby Isaac   Not collective
8959be51f97SToby Isaac 
8969be51f97SToby Isaac   Input Parameter:
8979be51f97SToby Isaac . dm - the forest
8989be51f97SToby Isaac 
8999be51f97SToby Isaac   Output Parameter:
9009be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL)
9019be51f97SToby Isaac 
9029be51f97SToby Isaac   Level: advanced
9039be51f97SToby Isaac 
9049be51f97SToby Isaac .seealso: DMForestSetAdaptivityStrategy()
9059be51f97SToby Isaac @*/
906c7eeac06SToby Isaac PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
907c7eeac06SToby Isaac {
908c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
909c7eeac06SToby Isaac 
910c7eeac06SToby Isaac   PetscFunctionBegin;
911c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
912c7eeac06SToby Isaac   PetscValidPointer(adaptStrategy,2);
913c7eeac06SToby Isaac   *adaptStrategy = forest->adaptStrategy;
914c7eeac06SToby Isaac   PetscFunctionReturn(0);
915c7eeac06SToby Isaac }
916c7eeac06SToby Isaac 
9172a133e43SToby Isaac /*@
9182a133e43SToby Isaac   DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning,
9192a133e43SToby Isaac   etc.) was successful.  PETSC_FALSE indicates that the post-adaptation forest is the same as the pre-adpatation
9202a133e43SToby Isaac   forest.  A requested adaptation may have been unsuccessful if, for example, the requested refinement would have
9212a133e43SToby Isaac   exceeded the maximum refinement level.
9222a133e43SToby Isaac 
9232a133e43SToby Isaac   Collective on dm
9242a133e43SToby Isaac 
9252a133e43SToby Isaac   Input Parameter:
9262a133e43SToby Isaac 
9272a133e43SToby Isaac . dm - the post-adaptation forest
9282a133e43SToby Isaac 
9292a133e43SToby Isaac   Output Parameter:
9302a133e43SToby Isaac 
9312a133e43SToby Isaac . success - PETSC_TRUE if the post-adaptation forest is different from the pre-adaptation forest.
9322a133e43SToby Isaac 
9332a133e43SToby Isaac   Level: intermediate
9342a133e43SToby Isaac 
9352a133e43SToby Isaac .see
9362a133e43SToby Isaac @*/
9372a133e43SToby Isaac PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success)
9382a133e43SToby Isaac {
9392a133e43SToby Isaac   DM_Forest      *forest;
9402a133e43SToby Isaac   PetscErrorCode ierr;
9412a133e43SToby Isaac 
9422a133e43SToby Isaac   PetscFunctionBegin;
9432a133e43SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9442a133e43SToby Isaac   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() has not been called yet.");
9452a133e43SToby Isaac   forest = (DM_Forest *) dm->data;
9462a133e43SToby Isaac   ierr = (forest->getadaptivitysuccess)(dm,success);CHKERRQ(ierr);
9472a133e43SToby Isaac   PetscFunctionReturn(0);
9482a133e43SToby Isaac }
9492a133e43SToby Isaac 
950bf9b5d84SToby Isaac /*@
951bf9b5d84SToby Isaac   DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed
952bf9b5d84SToby 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().
953bf9b5d84SToby Isaac 
954bf9b5d84SToby Isaac   Logically collective on dm
955bf9b5d84SToby Isaac 
956bf9b5d84SToby Isaac   Input Parameters:
957bf9b5d84SToby Isaac + dm - the post-adaptation forest
958bf9b5d84SToby Isaac - computeSF - default PETSC_TRUE
959bf9b5d84SToby Isaac 
960bf9b5d84SToby Isaac   Level: advanced
961bf9b5d84SToby Isaac 
962bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
963bf9b5d84SToby Isaac @*/
964bf9b5d84SToby Isaac PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
965bf9b5d84SToby Isaac {
966bf9b5d84SToby Isaac   DM_Forest *forest;
967bf9b5d84SToby Isaac 
968bf9b5d84SToby Isaac   PetscFunctionBegin;
969bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
970bf9b5d84SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called");
971bf9b5d84SToby Isaac   forest                 = (DM_Forest*) dm->data;
972bf9b5d84SToby Isaac   forest->computeAdaptSF = computeSF;
973bf9b5d84SToby Isaac   PetscFunctionReturn(0);
974bf9b5d84SToby Isaac }
975bf9b5d84SToby Isaac 
9760eb7e1eaSToby Isaac PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time)
97780b27e07SToby Isaac {
97880b27e07SToby Isaac   DM_Forest      *forest;
97980b27e07SToby Isaac   PetscErrorCode ierr;
98080b27e07SToby Isaac 
98180b27e07SToby Isaac   PetscFunctionBegin;
98280b27e07SToby Isaac   PetscValidHeaderSpecific(dmIn   ,DM_CLASSID  ,1);
98380b27e07SToby Isaac   PetscValidHeaderSpecific(vecIn  ,VEC_CLASSID ,2);
98480b27e07SToby Isaac   PetscValidHeaderSpecific(dmOut  ,DM_CLASSID  ,3);
98580b27e07SToby Isaac   PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,4);
98680b27e07SToby Isaac   forest = (DM_Forest *) dmIn->data;
98780b27e07SToby Isaac   if (!forest->transfervec) SETERRQ(PetscObjectComm((PetscObject)dmIn),PETSC_ERR_SUP,"DMForestTransferVec() not implemented");
9880eb7e1eaSToby Isaac   ierr = (forest->transfervec)(dmIn,vecIn,dmOut,vecOut,useBCs,time);CHKERRQ(ierr);
98980b27e07SToby Isaac   PetscFunctionReturn(0);
99080b27e07SToby Isaac }
99180b27e07SToby Isaac 
992bf9b5d84SToby Isaac /*@
993bf9b5d84SToby Isaac   DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the
994bf9b5d84SToby Isaac   pre-adaptation forest to the post-adaptiation forest.  After DMSetUp() is called, these transfer PetscSFs can be
995bf9b5d84SToby Isaac   accessed with DMForestGetAdaptivitySF().
996bf9b5d84SToby Isaac 
997bf9b5d84SToby Isaac   Not collective
998bf9b5d84SToby Isaac 
999bf9b5d84SToby Isaac   Input Parameter:
1000bf9b5d84SToby Isaac . dm - the post-adaptation forest
1001bf9b5d84SToby Isaac 
1002bf9b5d84SToby Isaac   Output Parameter:
1003bf9b5d84SToby Isaac . computeSF - default PETSC_TRUE
1004bf9b5d84SToby Isaac 
1005bf9b5d84SToby Isaac   Level: advanced
1006bf9b5d84SToby Isaac 
1007bf9b5d84SToby Isaac .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
1008bf9b5d84SToby Isaac @*/
1009bf9b5d84SToby Isaac PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1010bf9b5d84SToby Isaac {
1011bf9b5d84SToby Isaac   DM_Forest *forest;
1012bf9b5d84SToby Isaac 
1013bf9b5d84SToby Isaac   PetscFunctionBegin;
1014bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1015bf9b5d84SToby Isaac   forest     = (DM_Forest*) dm->data;
1016bf9b5d84SToby Isaac   *computeSF = forest->computeAdaptSF;
1017bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1018bf9b5d84SToby Isaac }
1019bf9b5d84SToby Isaac 
1020bf9b5d84SToby Isaac /*@
1021bf9b5d84SToby Isaac   DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest.
1022bf9b5d84SToby Isaac   Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be
1023bf9b5d84SToby Isaac   some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa.  Therefore there are two
1024bf9b5d84SToby Isaac   PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates
1025bf9b5d84SToby Isaac   pre-adaptation fine cells to post-adaptation coarse cells.
1026bf9b5d84SToby Isaac 
1027bf9b5d84SToby Isaac   Not collective
1028bf9b5d84SToby Isaac 
1029bf9b5d84SToby Isaac   Input Parameter:
1030bf9b5d84SToby Isaac   dm - the post-adaptation forest
1031bf9b5d84SToby Isaac 
1032bf9b5d84SToby Isaac   Output Parameter:
10330f17b9e3SToby Isaac   preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
10340f17b9e3SToby Isaac   coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-
1035bf9b5d84SToby Isaac 
1036bf9b5d84SToby Isaac   Level: advanced
1037bf9b5d84SToby Isaac 
1038bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF()
1039bf9b5d84SToby Isaac @*/
10400f17b9e3SToby Isaac PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1041bf9b5d84SToby Isaac {
1042bf9b5d84SToby Isaac   DM_Forest      *forest;
1043bf9b5d84SToby Isaac   PetscErrorCode ierr;
1044bf9b5d84SToby Isaac 
1045bf9b5d84SToby Isaac   PetscFunctionBegin;
1046bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1047bf9b5d84SToby Isaac   ierr   = DMSetUp(dm);CHKERRQ(ierr);
1048bf9b5d84SToby Isaac   forest = (DM_Forest*) dm->data;
1049f885a11aSToby Isaac   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1050f885a11aSToby Isaac   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
1051bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1052bf9b5d84SToby Isaac }
1053bf9b5d84SToby Isaac 
10549be51f97SToby Isaac /*@
10559be51f97SToby Isaac   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
10569be51f97SToby Isaac   indicate that the diameter of neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may
10579be51f97SToby Isaac   only support one particular choice of grading factor.
10589be51f97SToby Isaac 
10599be51f97SToby Isaac   Logically collective on dm
10609be51f97SToby Isaac 
10619be51f97SToby Isaac   Input Parameters:
10629be51f97SToby Isaac + dm - the forest
10639be51f97SToby Isaac - grade - the grading factor
10649be51f97SToby Isaac 
10659be51f97SToby Isaac   Level: advanced
10669be51f97SToby Isaac 
10679be51f97SToby Isaac .seealso: DMForestGetGradeFactor()
10689be51f97SToby Isaac @*/
1069c7eeac06SToby Isaac PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1070c7eeac06SToby Isaac {
1071c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1072c7eeac06SToby Isaac 
1073c7eeac06SToby Isaac   PetscFunctionBegin;
1074c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1075ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup");
1076c7eeac06SToby Isaac   forest->gradeFactor = grade;
1077c7eeac06SToby Isaac   PetscFunctionReturn(0);
1078c7eeac06SToby Isaac }
1079c7eeac06SToby Isaac 
10809be51f97SToby Isaac /*@
10819be51f97SToby Isaac   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
10829be51f97SToby Isaac   neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may only support one particular
10839be51f97SToby Isaac   choice of grading factor.
10849be51f97SToby Isaac 
10859be51f97SToby Isaac   Not collective
10869be51f97SToby Isaac 
10879be51f97SToby Isaac   Input Parameter:
10889be51f97SToby Isaac . dm - the forest
10899be51f97SToby Isaac 
10909be51f97SToby Isaac   Output Parameter:
10919be51f97SToby Isaac . grade - the grading factor
10929be51f97SToby Isaac 
10939be51f97SToby Isaac   Level: advanced
10949be51f97SToby Isaac 
10959be51f97SToby Isaac .seealso: DMForestSetGradeFactor()
10969be51f97SToby Isaac @*/
1097c7eeac06SToby Isaac PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1098c7eeac06SToby Isaac {
1099c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1100c7eeac06SToby Isaac 
1101c7eeac06SToby Isaac   PetscFunctionBegin;
1102c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1103c7eeac06SToby Isaac   PetscValidIntPointer(grade,2);
1104c7eeac06SToby Isaac   *grade = forest->gradeFactor;
1105c7eeac06SToby Isaac   PetscFunctionReturn(0);
1106c7eeac06SToby Isaac }
1107c7eeac06SToby Isaac 
11089be51f97SToby Isaac /*@
11099be51f97SToby Isaac   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
11109be51f97SToby Isaac   the cell weight (see DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be
11119be51f97SToby Isaac   (cellWeight) * (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on
11129be51f97SToby Isaac   its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
11139be51f97SToby Isaac   computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11149be51f97SToby Isaac 
11159be51f97SToby Isaac   Logically collective on dm
11169be51f97SToby Isaac 
11179be51f97SToby Isaac   Input Parameters:
11189be51f97SToby Isaac + dm - the forest
11199be51f97SToby Isaac - weightsFactors - default 1.
11209be51f97SToby Isaac 
11219be51f97SToby Isaac   Level: advanced
11229be51f97SToby Isaac 
11239be51f97SToby Isaac .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights()
11249be51f97SToby Isaac @*/
1125ef51cf95SToby Isaac PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
1126c7eeac06SToby Isaac {
1127c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1128c7eeac06SToby Isaac 
1129c7eeac06SToby Isaac   PetscFunctionBegin;
1130c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1131ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup");
1132c7eeac06SToby Isaac   forest->weightsFactor = weightsFactor;
1133c7eeac06SToby Isaac   PetscFunctionReturn(0);
1134c7eeac06SToby Isaac }
1135c7eeac06SToby Isaac 
11369be51f97SToby Isaac /*@
11379be51f97SToby Isaac   DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
11389be51f97SToby Isaac   DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be (cellWeight) *
11399be51f97SToby Isaac   (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on its level; a
11409be51f97SToby Isaac   factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation
11419be51f97SToby Isaac   associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11429be51f97SToby Isaac 
11439be51f97SToby Isaac   Not collective
11449be51f97SToby Isaac 
11459be51f97SToby Isaac   Input Parameter:
11469be51f97SToby Isaac . dm - the forest
11479be51f97SToby Isaac 
11489be51f97SToby Isaac   Output Parameter:
11499be51f97SToby Isaac . weightsFactors - default 1.
11509be51f97SToby Isaac 
11519be51f97SToby Isaac   Level: advanced
11529be51f97SToby Isaac 
11539be51f97SToby Isaac .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights()
11549be51f97SToby Isaac @*/
1155ef51cf95SToby Isaac PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1156c7eeac06SToby Isaac {
1157c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1158c7eeac06SToby Isaac 
1159c7eeac06SToby Isaac   PetscFunctionBegin;
1160c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1161c7eeac06SToby Isaac   PetscValidRealPointer(weightsFactor,2);
1162c7eeac06SToby Isaac   *weightsFactor = forest->weightsFactor;
1163c7eeac06SToby Isaac   PetscFunctionReturn(0);
1164c7eeac06SToby Isaac }
1165c7eeac06SToby Isaac 
11669be51f97SToby Isaac /*@
11679be51f97SToby Isaac   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
11689be51f97SToby Isaac 
11699be51f97SToby Isaac   Not collective
11709be51f97SToby Isaac 
11719be51f97SToby Isaac   Input Parameter:
11729be51f97SToby Isaac . dm - the forest
11739be51f97SToby Isaac 
11749be51f97SToby Isaac   Output Parameters:
11759be51f97SToby Isaac + cStart - the first cell on this process
11769be51f97SToby Isaac - cEnd - one after the final cell on this process
11779be51f97SToby Isaac 
11781a244344SSatish Balay   Level: intermediate
11799be51f97SToby Isaac 
11809be51f97SToby Isaac .seealso: DMForestGetCellSF()
11819be51f97SToby Isaac @*/
1182c7eeac06SToby Isaac PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1183c7eeac06SToby Isaac {
1184c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1185c7eeac06SToby Isaac   PetscErrorCode ierr;
1186c7eeac06SToby Isaac 
1187c7eeac06SToby Isaac   PetscFunctionBegin;
1188c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1189c7eeac06SToby Isaac   PetscValidIntPointer(cStart,2);
1190c7eeac06SToby Isaac   PetscValidIntPointer(cEnd,2);
1191c7eeac06SToby Isaac   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) {
1192c7eeac06SToby Isaac     ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr);
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 
12129be51f97SToby Isaac .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   PetscErrorCode ierr;
1218c7eeac06SToby Isaac 
1219c7eeac06SToby Isaac   PetscFunctionBegin;
1220c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1221c7eeac06SToby Isaac   PetscValidPointer(cellSF,2);
1222c7eeac06SToby Isaac   if ((!forest->cellSF) && forest->createcellsf) {
1223c7eeac06SToby Isaac     ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr);
1224c7eeac06SToby Isaac   }
1225c7eeac06SToby Isaac   *cellSF = forest->cellSF;
1226c7eeac06SToby Isaac   PetscFunctionReturn(0);
1227c7eeac06SToby Isaac }
1228c7eeac06SToby Isaac 
12299be51f97SToby Isaac /*@C
12309be51f97SToby Isaac   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
12319be51f97SToby Isaac   DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination).  The
1232cd3c525cSToby Isaac   interpretation of the label values is up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP,
1233cd3c525cSToby Isaac   DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have been reserved as choices that should be accepted by all subtypes.
12349be51f97SToby Isaac 
12359be51f97SToby Isaac   Logically collective on dm
12369be51f97SToby Isaac 
12379be51f97SToby Isaac   Input Parameters:
12389be51f97SToby Isaac - dm - the forest
1239a1b0c543SToby Isaac + adaptLabel - the label in the pre-adaptation forest
12409be51f97SToby Isaac 
12419be51f97SToby Isaac   Level: intermediate
12429be51f97SToby Isaac 
12439be51f97SToby Isaac .seealso DMForestGetAdaptivityLabel()
12449be51f97SToby Isaac @*/
1245a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel)
1246c7eeac06SToby Isaac {
1247c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1248c7eeac06SToby Isaac   PetscErrorCode ierr;
1249c7eeac06SToby Isaac 
1250c7eeac06SToby Isaac   PetscFunctionBegin;
1251c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1252a1b0c543SToby Isaac   adaptLabel->refct++;
1253a1b0c543SToby Isaac   if (forest->adaptLabel) {ierr = DMLabelDestroy(&forest->adaptLabel);CHKERRQ(ierr);}
1254a1b0c543SToby Isaac   forest->adaptLabel = adaptLabel;
1255c7eeac06SToby Isaac   PetscFunctionReturn(0);
1256c7eeac06SToby Isaac }
1257c7eeac06SToby Isaac 
12589be51f97SToby Isaac /*@C
12599be51f97SToby Isaac   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that
12609be51f97SToby Isaac   holds the adaptation flags (refinement, coarsening, or some combination).  The interpretation of the label values is
1261cd3c525cSToby Isaac   up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have
1262cd3c525cSToby Isaac   been reserved as choices that should be accepted by all subtypes.
12639be51f97SToby Isaac 
12649be51f97SToby Isaac   Not collective
12659be51f97SToby Isaac 
12669be51f97SToby Isaac   Input Parameter:
12679be51f97SToby Isaac . dm - the forest
12689be51f97SToby Isaac 
12699be51f97SToby Isaac   Output Parameter:
12709be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest
12719be51f97SToby Isaac 
12729be51f97SToby Isaac   Level: intermediate
12739be51f97SToby Isaac 
12749be51f97SToby Isaac .seealso DMForestSetAdaptivityLabel()
12759be51f97SToby Isaac @*/
1276a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel)
1277c7eeac06SToby Isaac {
1278c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1279c7eeac06SToby Isaac 
1280c7eeac06SToby Isaac   PetscFunctionBegin;
1281c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1282ba936b91SToby Isaac   *adaptLabel = forest->adaptLabel;
1283c7eeac06SToby Isaac   PetscFunctionReturn(0);
1284c7eeac06SToby Isaac }
1285c7eeac06SToby Isaac 
12869be51f97SToby Isaac /*@
12879be51f97SToby Isaac   DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
12889be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
12899be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
12909be51f97SToby Isaac   of 1.
12919be51f97SToby Isaac 
12929be51f97SToby Isaac   Logically collective on dm
12939be51f97SToby Isaac 
12949be51f97SToby Isaac   Input Parameters:
12959be51f97SToby Isaac + dm - the forest
12969be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
12979be51f97SToby Isaac - copyMode - how weights should reference weights
12989be51f97SToby Isaac 
12999be51f97SToby Isaac   Level: advanced
13009be51f97SToby Isaac 
13019be51f97SToby Isaac .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity()
13029be51f97SToby Isaac @*/
1303c7eeac06SToby Isaac PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1304c7eeac06SToby Isaac {
1305c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1306c7eeac06SToby Isaac   PetscInt       cStart, cEnd;
1307c7eeac06SToby Isaac   PetscErrorCode ierr;
1308c7eeac06SToby Isaac 
1309c7eeac06SToby Isaac   PetscFunctionBegin;
1310c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1311c7eeac06SToby Isaac   ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr);
1312c7eeac06SToby Isaac   if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd);
1313c7eeac06SToby Isaac   if (copyMode == PETSC_COPY_VALUES) {
1314c7eeac06SToby Isaac     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) {
1315c7eeac06SToby Isaac       ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr);
1316c7eeac06SToby Isaac     }
1317c7eeac06SToby Isaac     ierr                        = PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));CHKERRQ(ierr);
1318c7eeac06SToby Isaac     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
1319c7eeac06SToby Isaac     PetscFunctionReturn(0);
1320c7eeac06SToby Isaac   }
1321c7eeac06SToby Isaac   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) {
1322c7eeac06SToby Isaac     ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr);
1323c7eeac06SToby Isaac   }
1324c7eeac06SToby Isaac   forest->cellWeights         = weights;
1325c7eeac06SToby Isaac   forest->cellWeightsCopyMode = copyMode;
1326c7eeac06SToby Isaac   PetscFunctionReturn(0);
1327c7eeac06SToby Isaac }
1328c7eeac06SToby Isaac 
13299be51f97SToby Isaac /*@
13309be51f97SToby Isaac   DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
13319be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
13329be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
13339be51f97SToby Isaac   of 1.
13349be51f97SToby Isaac 
13359be51f97SToby Isaac   Not collective
13369be51f97SToby Isaac 
13379be51f97SToby Isaac   Input Parameter:
13389be51f97SToby Isaac . dm - the forest
13399be51f97SToby Isaac 
13409be51f97SToby Isaac   Output Parameter:
13419be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
13429be51f97SToby Isaac 
13439be51f97SToby Isaac   Level: advanced
13449be51f97SToby Isaac 
13459be51f97SToby Isaac .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity()
13469be51f97SToby Isaac @*/
1347c7eeac06SToby Isaac PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1348c7eeac06SToby Isaac {
1349c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1350c7eeac06SToby Isaac 
1351c7eeac06SToby Isaac   PetscFunctionBegin;
1352c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1353c7eeac06SToby Isaac   PetscValidPointer(weights,2);
1354c7eeac06SToby Isaac   *weights = forest->cellWeights;
1355c7eeac06SToby Isaac   PetscFunctionReturn(0);
1356c7eeac06SToby Isaac }
1357c7eeac06SToby Isaac 
13589be51f97SToby Isaac /*@
13599be51f97SToby Isaac   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
13609be51f97SToby Isaac   a pre-adaptation forest (see DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each
13619be51f97SToby Isaac   process's cells to the process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that
13629be51f97SToby Isaac   the current process should not have any cells after repartitioning.
13639be51f97SToby Isaac 
13649be51f97SToby Isaac   Logically Collective on dm
13659be51f97SToby Isaac 
13669be51f97SToby Isaac   Input parameters:
13679be51f97SToby Isaac + dm - the forest
13689be51f97SToby Isaac - capacity - this process's capacity
13699be51f97SToby Isaac 
13709be51f97SToby Isaac   Level: advanced
13719be51f97SToby Isaac 
13729be51f97SToby Isaac .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
13739be51f97SToby Isaac @*/
1374c7eeac06SToby Isaac PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1375c7eeac06SToby Isaac {
1376c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1377c7eeac06SToby Isaac 
1378c7eeac06SToby Isaac   PetscFunctionBegin;
1379c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1380ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup");
1381c7eeac06SToby Isaac   if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity);
1382c7eeac06SToby Isaac   forest->weightCapacity = capacity;
1383c7eeac06SToby Isaac   PetscFunctionReturn(0);
1384c7eeac06SToby Isaac }
1385c7eeac06SToby Isaac 
13869be51f97SToby Isaac /*@
13879be51f97SToby Isaac   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
13889be51f97SToby Isaac   DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each process's cells to the
13899be51f97SToby Isaac   process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that the current process
13909be51f97SToby Isaac   should not have any cells after repartitioning.
13919be51f97SToby Isaac 
13929be51f97SToby Isaac   Not collective
13939be51f97SToby Isaac 
13949be51f97SToby Isaac   Input parameter:
13959be51f97SToby Isaac . dm - the forest
13969be51f97SToby Isaac 
13979be51f97SToby Isaac   Output parameter:
13989be51f97SToby Isaac . capacity - this process's capacity
13999be51f97SToby Isaac 
14009be51f97SToby Isaac   Level: advanced
14019be51f97SToby Isaac 
14029be51f97SToby Isaac .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
14039be51f97SToby Isaac @*/
1404c7eeac06SToby Isaac PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1405c7eeac06SToby Isaac {
1406c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1407c7eeac06SToby Isaac 
1408c7eeac06SToby Isaac   PetscFunctionBegin;
1409c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1410c7eeac06SToby Isaac   PetscValidRealPointer(capacity,2);
1411c7eeac06SToby Isaac   *capacity = forest->weightCapacity;
1412c7eeac06SToby Isaac   PetscFunctionReturn(0);
1413c7eeac06SToby Isaac }
1414c7eeac06SToby Isaac 
14150709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm)
1416db4d5e8cSToby Isaac {
1417db4d5e8cSToby Isaac   DM_Forest                  *forest = (DM_Forest*) dm->data;
141856ba9f64SToby Isaac   PetscBool                  flg, flg1, flg2, flg3, flg4;
1419dd8e54a2SToby Isaac   DMForestTopology           oldTopo;
1420c7eeac06SToby Isaac   char                       stringBuffer[256];
1421dd8e54a2SToby Isaac   PetscViewer                viewer;
1422dd8e54a2SToby Isaac   PetscViewerFormat          format;
142356ba9f64SToby Isaac   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1424c7eeac06SToby Isaac   PetscReal                  weightsFactor;
1425c7eeac06SToby Isaac   DMForestAdaptivityStrategy adaptStrategy;
1426db4d5e8cSToby Isaac   PetscErrorCode             ierr;
1427db4d5e8cSToby Isaac 
1428db4d5e8cSToby Isaac   PetscFunctionBegin;
1429db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
143058762b62SToby Isaac   forest->setfromoptionscalled = PETSC_TRUE;
1431dd8e54a2SToby Isaac   ierr                         = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr);
1432a3eda1eaSToby Isaac   ierr                         = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr);
143356ba9f64SToby Isaac   ierr                         = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);CHKERRQ(ierr);
143456ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr);
143556ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr);
143656ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr);
1437f885a11aSToby Isaac   if ((PetscInt) flg1 + (PetscInt) flg2 + (PetscInt) flg3 + (PetscInt) flg4 > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Specify only one of -dm_forest_{topology,base_dm,coarse_forest,fine_forest}");
143856ba9f64SToby Isaac   if (flg1) {
143956ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr);
144056ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
144120e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
144256ba9f64SToby Isaac   }
144356ba9f64SToby Isaac   if (flg2) {
1444dd8e54a2SToby Isaac     DM base;
1445dd8e54a2SToby Isaac 
1446dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr);
1447dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1448dd8e54a2SToby Isaac     ierr = DMLoad(base,viewer);CHKERRQ(ierr);
1449dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1450dd8e54a2SToby Isaac     ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr);
1451dd8e54a2SToby Isaac     ierr = DMDestroy(&base);CHKERRQ(ierr);
145256ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
145320e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
1454dd8e54a2SToby Isaac   }
145556ba9f64SToby Isaac   if (flg3) {
1456dd8e54a2SToby Isaac     DM coarse;
1457dd8e54a2SToby Isaac 
1458dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr);
1459dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1460dd8e54a2SToby Isaac     ierr = DMLoad(coarse,viewer);CHKERRQ(ierr);
1461dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
146220e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr);
1463dd8e54a2SToby Isaac     ierr = DMDestroy(&coarse);CHKERRQ(ierr);
146456ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
146556ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1466dd8e54a2SToby Isaac   }
146756ba9f64SToby Isaac   if (flg4) {
1468dd8e54a2SToby Isaac     DM fine;
1469dd8e54a2SToby Isaac 
1470dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr);
1471dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1472dd8e54a2SToby Isaac     ierr = DMLoad(fine,viewer);CHKERRQ(ierr);
1473dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
147420e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr);
1475dd8e54a2SToby Isaac     ierr = DMDestroy(&fine);CHKERRQ(ierr);
147656ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
147756ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1478dd8e54a2SToby Isaac   }
1479dd8e54a2SToby Isaac   ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr);
1480dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);CHKERRQ(ierr);
1481dd8e54a2SToby Isaac   if (flg) {
1482dd8e54a2SToby Isaac     ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr);
1483f885a11aSToby Isaac   } else {
1484dd8e54a2SToby Isaac     ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr);
1485dd8e54a2SToby Isaac     ierr = PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);CHKERRQ(ierr);
1486dd8e54a2SToby Isaac     if (flg) {
1487dd8e54a2SToby Isaac       ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr);
1488dd8e54a2SToby Isaac     }
1489dd8e54a2SToby Isaac   }
1490dd8e54a2SToby Isaac   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
1491dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);CHKERRQ(ierr);
1492dd8e54a2SToby Isaac   if (flg) {
1493dd8e54a2SToby Isaac     ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr);
1494dd8e54a2SToby Isaac   }
1495a6121fbdSMatthew G. Knepley #if 0
1496a6121fbdSMatthew G. Knepley   ierr = PetscOptionsInt("-dm_refine","equivalent to -dm_forest_set_minimum_refinement and -dm_forest_set_initial_refinement with the same value",NULL,minRefinement,&minRefinement,&flg);CHKERRQ(ierr);
1497a6121fbdSMatthew G. Knepley   if (flg) {
1498a6121fbdSMatthew G. Knepley     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1499a6121fbdSMatthew G. Knepley     ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr);
1500a6121fbdSMatthew G. Knepley   }
1501a6121fbdSMatthew G. Knepley   ierr = PetscOptionsInt("-dm_refine_hierarchy","equivalent to -dm_forest_set_minimum_refinement 0 and -dm_forest_set_initial_refinement",NULL,initRefinement,&initRefinement,&flg);CHKERRQ(ierr);
1502a6121fbdSMatthew G. Knepley   if (flg) {
1503a6121fbdSMatthew G. Knepley     ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr);
1504a6121fbdSMatthew G. Knepley     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
1505a6121fbdSMatthew G. Knepley   }
1506a6121fbdSMatthew G. Knepley #endif
1507dd8e54a2SToby Isaac   ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr);
1508dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);CHKERRQ(ierr);
1509dd8e54a2SToby Isaac   if (flg) {
1510dd8e54a2SToby Isaac     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1511db4d5e8cSToby Isaac   }
151256ba9f64SToby Isaac   ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr);
151356ba9f64SToby Isaac   ierr = PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);CHKERRQ(ierr);
151456ba9f64SToby Isaac   if (flg) {
151556ba9f64SToby Isaac     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
151656ba9f64SToby Isaac   }
1517c7eeac06SToby Isaac   ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr);
1518c7eeac06SToby Isaac   ierr = PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);CHKERRQ(ierr);
1519c7eeac06SToby Isaac   if (flg) {
1520c7eeac06SToby Isaac     ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr);
1521c7eeac06SToby Isaac   }
1522c7eeac06SToby Isaac   ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr);
1523c7eeac06SToby Isaac   ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);CHKERRQ(ierr);
1524c7eeac06SToby Isaac   if (flg) {
1525c7eeac06SToby Isaac     ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr);
1526c7eeac06SToby Isaac   }
1527c7eeac06SToby Isaac   ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr);
1528c7eeac06SToby Isaac   ierr = PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);CHKERRQ(ierr);
1529c7eeac06SToby Isaac   if (flg) {
1530c7eeac06SToby Isaac     ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr);
1531c7eeac06SToby Isaac   }
1532c7eeac06SToby Isaac   ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr);
1533c7eeac06SToby Isaac   ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr);
1534c7eeac06SToby Isaac   if (flg) {
1535c7eeac06SToby Isaac     ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr);
1536c7eeac06SToby Isaac   }
1537db4d5e8cSToby Isaac   ierr = PetscOptionsTail();CHKERRQ(ierr);
1538db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1539db4d5e8cSToby Isaac }
1540db4d5e8cSToby Isaac 
1541d8984e3bSMatthew G. Knepley PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
1542d8984e3bSMatthew G. Knepley {
1543d8984e3bSMatthew G. Knepley   PetscErrorCode ierr;
1544d8984e3bSMatthew G. Knepley 
1545d8984e3bSMatthew G. Knepley   PetscFunctionBegin;
1546d8984e3bSMatthew G. Knepley   if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);}
1547d8984e3bSMatthew G. Knepley   ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
1548d8984e3bSMatthew G. Knepley   PetscFunctionReturn(0);
1549d8984e3bSMatthew G. Knepley }
1550d8984e3bSMatthew G. Knepley 
15515421bac9SToby Isaac PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
15525421bac9SToby Isaac {
15535421bac9SToby Isaac   DMLabel        refine;
15545421bac9SToby Isaac   DM             fineDM;
15555421bac9SToby Isaac   PetscErrorCode ierr;
15565421bac9SToby Isaac 
15575421bac9SToby Isaac   PetscFunctionBegin;
15585421bac9SToby Isaac   ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr);
15595421bac9SToby Isaac   if (fineDM) {
15605421bac9SToby Isaac     ierr       = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr);
15615421bac9SToby Isaac     *dmRefined = fineDM;
15625421bac9SToby Isaac     PetscFunctionReturn(0);
15635421bac9SToby Isaac   }
15645421bac9SToby Isaac   ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr);
15655421bac9SToby Isaac   ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr);
15665421bac9SToby Isaac   if (!refine) {
1567a1b0c543SToby Isaac     ierr = DMLabelCreate("refine",&refine);CHKERRQ(ierr);
1568a1b0c543SToby Isaac     ierr = DMLabelSetDefaultValue(refine,DM_ADAPT_REFINE);CHKERRQ(ierr);
15695421bac9SToby Isaac   }
1570a1b0c543SToby Isaac   else {
1571a1b0c543SToby Isaac     refine->refct++;
1572a1b0c543SToby Isaac   }
1573a1b0c543SToby Isaac   ierr = DMForestSetAdaptivityLabel(*dmRefined,refine);CHKERRQ(ierr);
1574a1b0c543SToby Isaac   ierr = DMLabelDestroy(&refine);CHKERRQ(ierr);
15755421bac9SToby Isaac   PetscFunctionReturn(0);
15765421bac9SToby Isaac }
15775421bac9SToby Isaac 
15785421bac9SToby Isaac PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
15795421bac9SToby Isaac {
15805421bac9SToby Isaac   DMLabel        coarsen;
15815421bac9SToby Isaac   DM             coarseDM;
15825421bac9SToby Isaac   PetscErrorCode ierr;
15835421bac9SToby Isaac 
15845421bac9SToby Isaac   PetscFunctionBegin;
15854098eed7SToby Isaac   {
15864098eed7SToby Isaac     PetscMPIInt mpiComparison;
15874098eed7SToby Isaac     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);
15884098eed7SToby Isaac 
15894098eed7SToby Isaac     ierr = MPI_Comm_compare(comm, dmcomm, &mpiComparison);CHKERRQ(ierr);
1590f885a11aSToby Isaac     if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet");
15914098eed7SToby Isaac   }
15925421bac9SToby Isaac   ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr);
15935421bac9SToby Isaac   if (coarseDM) {
15945421bac9SToby Isaac     ierr         = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr);
15955421bac9SToby Isaac     *dmCoarsened = coarseDM;
15965421bac9SToby Isaac     PetscFunctionReturn(0);
15975421bac9SToby Isaac   }
15985421bac9SToby Isaac   ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr);
1599a1b0c543SToby Isaac   ierr = DMForestSetAdaptivityPurpose(coarseDM,DM_ADAPT_COARSEN);CHKERRQ(ierr);
16005421bac9SToby Isaac   ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr);
16015421bac9SToby Isaac   if (!coarsen) {
1602a1b0c543SToby Isaac     ierr = DMLabelCreate("coarsen",&coarsen);CHKERRQ(ierr);
1603a1b0c543SToby Isaac     ierr = DMLabelSetDefaultValue(coarsen,DM_ADAPT_COARSEN);CHKERRQ(ierr);
1604a1b0c543SToby Isaac   } else {
1605a1b0c543SToby Isaac     coarsen->refct++;
16065421bac9SToby Isaac   }
1607a1b0c543SToby Isaac   ierr = DMForestSetAdaptivityLabel(*dmCoarsened,coarsen);CHKERRQ(ierr);
1608a1b0c543SToby Isaac   ierr = DMLabelDestroy(&coarsen);CHKERRQ(ierr);
16095421bac9SToby Isaac   PetscFunctionReturn(0);
16105421bac9SToby Isaac }
16115421bac9SToby Isaac 
1612a1b0c543SToby Isaac static PetscErrorCode DMAdaptLabel_Forest(DM dm, DMLabel label, DM *adaptedDM)
161309350103SToby Isaac {
161409350103SToby Isaac   PetscBool      success;
161509350103SToby Isaac   PetscErrorCode ierr;
161609350103SToby Isaac 
161709350103SToby Isaac   PetscFunctionBegin;
161809350103SToby Isaac   ierr = DMForestTemplate(dm,PetscObjectComm((PetscObject)dm),adaptedDM);CHKERRQ(ierr);
1619a1b0c543SToby Isaac   ierr = DMForestSetAdaptivityLabel(*adaptedDM,label);CHKERRQ(ierr);
162009350103SToby Isaac   ierr = DMSetUp(*adaptedDM);CHKERRQ(ierr);
162109350103SToby Isaac   ierr = DMForestGetAdaptivitySuccess(*adaptedDM,&success);CHKERRQ(ierr);
162209350103SToby Isaac   if (!success) {
162309350103SToby Isaac     ierr = DMDestroy(adaptedDM);CHKERRQ(ierr);
162409350103SToby Isaac     *adaptedDM = NULL;
162509350103SToby Isaac   }
162609350103SToby Isaac   PetscFunctionReturn(0);
162709350103SToby Isaac }
162809350103SToby Isaac 
1629d222f98bSToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm)
1630d222f98bSToby Isaac {
1631d222f98bSToby Isaac   PetscErrorCode ierr;
1632d222f98bSToby Isaac 
1633d222f98bSToby Isaac   PetscFunctionBegin;
1634d222f98bSToby Isaac   ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr);
1635d222f98bSToby Isaac 
1636d222f98bSToby Isaac   dm->ops->clone          = DMClone_Forest;
1637d222f98bSToby Isaac   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1638d222f98bSToby Isaac   dm->ops->destroy        = DMDestroy_Forest;
1639d8984e3bSMatthew G. Knepley   dm->ops->createsubdm    = DMCreateSubDM_Forest;
16405421bac9SToby Isaac   dm->ops->refine         = DMRefine_Forest;
16415421bac9SToby Isaac   dm->ops->coarsen        = DMCoarsen_Forest;
16420d1cd5e0SMatthew G. Knepley   dm->ops->adaptlabel     = DMAdaptLabel_Forest;
1643d222f98bSToby Isaac   PetscFunctionReturn(0);
1644d222f98bSToby Isaac }
1645d222f98bSToby Isaac 
16469be51f97SToby Isaac /*MC
16479be51f97SToby Isaac 
1648bae1f979SBarry Smith      DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh.  Forests usually have a base DM
1649bae1f979SBarry Smith   (see DMForestGetBaseDM()), from which it is refined.  The refinement and partitioning of forests is considered
1650bae1f979SBarry Smith   immutable after DMSetUp() is called.  To adapt a mesh, one should call DMForestTemplate() to create a new mesh that
1651bae1f979SBarry Smith   will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new
1652bae1f979SBarry Smith   mesh.
1653bae1f979SBarry Smith 
1654bae1f979SBarry Smith   To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the
1655bae1f979SBarry Smith   previous mesh whose values indicate which cells should be refined (DM_ADAPT_REFINE) or coarsened (DM_ADAPT_COARSEN)
1656bae1f979SBarry Smith   and how (subtypes are free to allow additional values for things like anisotropic refinement).  The label should be
1657bae1f979SBarry Smith   given to the *new* mesh with DMForestSetAdaptivityLabel().
16589be51f97SToby Isaac 
16599be51f97SToby Isaac   Level: advanced
16609be51f97SToby Isaac 
16619be51f97SToby Isaac .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel()
16629be51f97SToby Isaac M*/
16639be51f97SToby Isaac 
1664db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1665db4d5e8cSToby Isaac {
1666db4d5e8cSToby Isaac   DM_Forest      *forest;
1667db4d5e8cSToby Isaac   PetscErrorCode ierr;
1668db4d5e8cSToby Isaac 
1669db4d5e8cSToby Isaac   PetscFunctionBegin;
1670db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1671db4d5e8cSToby Isaac   ierr                         = PetscNewLog(dm,&forest);CHKERRQ(ierr);
1672db4d5e8cSToby Isaac   dm->dim                      = 0;
1673db4d5e8cSToby Isaac   dm->data                     = forest;
1674db4d5e8cSToby Isaac   forest->refct                = 1;
1675db4d5e8cSToby Isaac   forest->data                 = NULL;
167658762b62SToby Isaac   forest->setfromoptionscalled = PETSC_FALSE;
1677db4d5e8cSToby Isaac   forest->topology             = NULL;
1678cd3c525cSToby Isaac   forest->adapt                = NULL;
1679db4d5e8cSToby Isaac   forest->base                 = NULL;
16806a87ffbfSToby Isaac   forest->adaptPurpose         = DM_ADAPT_DETERMINE;
1681db4d5e8cSToby Isaac   forest->adjDim               = PETSC_DEFAULT;
1682db4d5e8cSToby Isaac   forest->overlap              = PETSC_DEFAULT;
1683db4d5e8cSToby Isaac   forest->minRefinement        = PETSC_DEFAULT;
1684db4d5e8cSToby Isaac   forest->maxRefinement        = PETSC_DEFAULT;
168556ba9f64SToby Isaac   forest->initRefinement       = PETSC_DEFAULT;
1686c7eeac06SToby Isaac   forest->cStart               = PETSC_DETERMINE;
1687c7eeac06SToby Isaac   forest->cEnd                 = PETSC_DETERMINE;
1688cd3c525cSToby Isaac   forest->cellSF               = NULL;
1689ebdf65a2SToby Isaac   forest->adaptLabel           = NULL;
1690db4d5e8cSToby Isaac   forest->gradeFactor          = 2;
1691db4d5e8cSToby Isaac   forest->cellWeights          = NULL;
1692db4d5e8cSToby Isaac   forest->cellWeightsCopyMode  = PETSC_USE_POINTER;
1693db4d5e8cSToby Isaac   forest->weightsFactor        = 1.;
1694db4d5e8cSToby Isaac   forest->weightCapacity       = 1.;
1695a73e2921SToby Isaac   ierr                         = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr);
1696d222f98bSToby Isaac   ierr                         = DMInitialize_Forest(dm);CHKERRQ(ierr);
1697db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1698db4d5e8cSToby Isaac }
1699db4d5e8cSToby Isaac 
1700