xref: /petsc/src/dm/impls/forest/forest.c (revision d8d19677bbccf95218448bee62e6b87f4513e133)
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;
13105e99e11SStefano Zampini   MatType                    mtype;
132a0452a8eSToby Isaac   PetscInt                   dim, overlap, ref, factor;
133a0452a8eSToby Isaac   DMForestAdaptivityStrategy strat;
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) {
1631fd50544SStefano Zampini     ierr = (*forest->ftemplate)(dm, *tdm);CHKERRQ(ierr);
164a0452a8eSToby Isaac   }
16520e8089bSToby Isaac   ierr = DMForestSetAdaptivityForest(*tdm,dm);CHKERRQ(ierr);
166e5e52638SMatthew G. Knepley   ierr = DMCopyDisc(dm,*tdm);CHKERRQ(ierr);
167795844e7SToby Isaac   ierr = DMGetApplicationContext(dm,&ctx);CHKERRQ(ierr);
168795844e7SToby Isaac   ierr = DMSetApplicationContext(*tdm,&ctx);CHKERRQ(ierr);
16990b157c4SStefano Zampini   {
17090b157c4SStefano Zampini     PetscBool            isper;
171795844e7SToby Isaac     const PetscReal      *maxCell, *L;
172795844e7SToby Isaac     const DMBoundaryType *bd;
173795844e7SToby Isaac 
17490b157c4SStefano Zampini     ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
17590b157c4SStefano Zampini     ierr = DMSetPeriodicity(*tdm,isper,maxCell,L,bd);CHKERRQ(ierr);
176795844e7SToby Isaac   }
17705e99e11SStefano Zampini   ierr = DMGetMatType(dm,&mtype);CHKERRQ(ierr);
17805e99e11SStefano Zampini   ierr = DMSetMatType(*tdm,mtype);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);
2061fd50544SStefano Zampini   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 
232fd292e60Sprj- .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 
28495452b02SPatrick Sanan   Notes:
28595452b02SPatrick Sanan     Currently the base DM must be a DMPLEX
286765b024eSBarry Smith 
2879be51f97SToby Isaac   Level: intermediate
2889be51f97SToby Isaac 
289fd292e60Sprj- .seealso: DMForestGetBaseDM()
2909be51f97SToby Isaac @*/
291dd8e54a2SToby Isaac PetscErrorCode DMForestSetBaseDM(DM dm, DM base)
292dd8e54a2SToby Isaac {
293dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
294dd8e54a2SToby Isaac   PetscInt       dim, dimEmbed;
295dd8e54a2SToby Isaac   PetscErrorCode ierr;
296dd8e54a2SToby Isaac 
297dd8e54a2SToby Isaac   PetscFunctionBegin;
298dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
299ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup");
300dd8e54a2SToby Isaac   ierr         = PetscObjectReference((PetscObject)base);CHKERRQ(ierr);
301dd8e54a2SToby Isaac   ierr         = DMDestroy(&forest->base);CHKERRQ(ierr);
302dd8e54a2SToby Isaac   forest->base = base;
303a0452a8eSToby Isaac   if (base) {
30428dfcf7cSStefano Zampini     PetscBool        isper;
30528dfcf7cSStefano Zampini     const PetscReal *maxCell, *L;
30628dfcf7cSStefano Zampini     const DMBoundaryType *bd;
30728dfcf7cSStefano Zampini 
308a0452a8eSToby Isaac     PetscValidHeaderSpecific(base, DM_CLASSID, 2);
309dd8e54a2SToby Isaac     ierr = DMGetDimension(base,&dim);CHKERRQ(ierr);
310dd8e54a2SToby Isaac     ierr = DMSetDimension(dm,dim);CHKERRQ(ierr);
311dd8e54a2SToby Isaac     ierr = DMGetCoordinateDim(base,&dimEmbed);CHKERRQ(ierr);
312dd8e54a2SToby Isaac     ierr = DMSetCoordinateDim(dm,dimEmbed);CHKERRQ(ierr);
31328dfcf7cSStefano Zampini     ierr = DMGetPeriodicity(base,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
31428dfcf7cSStefano Zampini     ierr = DMSetPeriodicity(dm,isper,maxCell,L,bd);CHKERRQ(ierr);
31528dfcf7cSStefano Zampini   } else {
31628dfcf7cSStefano Zampini     ierr = DMSetPeriodicity(dm,PETSC_FALSE,NULL,NULL,NULL);CHKERRQ(ierr);
317a0452a8eSToby Isaac   }
318dd8e54a2SToby Isaac   PetscFunctionReturn(0);
319dd8e54a2SToby Isaac }
320dd8e54a2SToby Isaac 
3219be51f97SToby Isaac /*@
3229be51f97SToby Isaac   DMForestGetBaseDM - Get the base DM of a DMForest forest.  The forest will be hierarchically refined from the base,
32368d54884SBarry Smith   and all refinements/coarsenings of the forest will share its base.  In general, two forest must share a base to be
3249be51f97SToby Isaac   comparable, to do things like construct interpolators.
3259be51f97SToby Isaac 
3269be51f97SToby Isaac   Not collective
3279be51f97SToby Isaac 
3289be51f97SToby Isaac   Input Parameter:
3299be51f97SToby Isaac . dm - the forest
3309be51f97SToby Isaac 
3319be51f97SToby Isaac   Output Parameter:
3329be51f97SToby Isaac . base - the base DM of the forest
3339be51f97SToby Isaac 
334367003a6SStefano Zampini   Notes:
335367003a6SStefano Zampini     After DMSetUp(), the base DM will be redundantly distributed across MPI processes
336367003a6SStefano Zampini 
3379be51f97SToby Isaac   Level: intermediate
3389be51f97SToby Isaac 
339fd292e60Sprj- .seealso: DMForestSetBaseDM()
3409be51f97SToby Isaac @*/
341dd8e54a2SToby Isaac PetscErrorCode DMForestGetBaseDM(DM dm, DM *base)
342dd8e54a2SToby Isaac {
343dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
344dd8e54a2SToby Isaac 
345dd8e54a2SToby Isaac   PetscFunctionBegin;
346dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
347dd8e54a2SToby Isaac   PetscValidPointer(base, 2);
348dd8e54a2SToby Isaac   *base = forest->base;
349dd8e54a2SToby Isaac   PetscFunctionReturn(0);
350dd8e54a2SToby Isaac }
351dd8e54a2SToby Isaac 
35299478f86SToby Isaac PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
353cf38a08cSToby Isaac {
354cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
355cf38a08cSToby Isaac 
356cf38a08cSToby Isaac   PetscFunctionBegin;
357cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
358cf38a08cSToby Isaac   forest->mapcoordinates    = func;
359cf38a08cSToby Isaac   forest->mapcoordinatesctx = ctx;
360cf38a08cSToby Isaac   PetscFunctionReturn(0);
361cf38a08cSToby Isaac }
362cf38a08cSToby Isaac 
36399478f86SToby Isaac PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func) (DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
364cf38a08cSToby Isaac {
365cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
366cf38a08cSToby Isaac 
367cf38a08cSToby Isaac   PetscFunctionBegin;
368cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
369cf38a08cSToby Isaac   if (func) *func = forest->mapcoordinates;
370cf38a08cSToby Isaac   if (ctx) *((void**) ctx) = forest->mapcoordinatesctx;
371cf38a08cSToby Isaac   PetscFunctionReturn(0);
372cf38a08cSToby Isaac }
373cf38a08cSToby Isaac 
3749be51f97SToby Isaac /*@
3759be51f97SToby Isaac   DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be
3769be51f97SToby Isaac   adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp().  Usually not needed
3779be51f97SToby Isaac   by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this
3789be51f97SToby Isaac   routine.
3799be51f97SToby Isaac 
380dffe73a3SToby Isaac   Note that this can be called after setup with adapt = NULL, which will clear all internal data related to the
381dffe73a3SToby Isaac   adaptivity forest from dm.  This way, repeatedly adapting does not leave stale DM objects in memory.
382dffe73a3SToby Isaac 
3839be51f97SToby Isaac   Logically collective on dm
3849be51f97SToby Isaac 
385*d8d19677SJose E. Roman   Input Parameters:
3869be51f97SToby Isaac + dm - the new forest, which will be constructed from adapt
3879be51f97SToby Isaac - adapt - the old forest
3889be51f97SToby Isaac 
3899be51f97SToby Isaac   Level: intermediate
3909be51f97SToby Isaac 
3919be51f97SToby Isaac .seealso: DMForestGetAdaptivityForest(), DMForestSetAdaptivityPurpose()
3929be51f97SToby Isaac @*/
393ba936b91SToby Isaac PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt)
394dd8e54a2SToby Isaac {
395dffe73a3SToby Isaac   DM_Forest      *forest, *adaptForest, *oldAdaptForest;
396dffe73a3SToby Isaac   DM             oldAdapt;
397456cc5b7SMatthew G. Knepley   PetscBool      isForest;
398dd8e54a2SToby Isaac   PetscErrorCode ierr;
399dd8e54a2SToby Isaac 
400dd8e54a2SToby Isaac   PetscFunctionBegin;
401dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4021fd50544SStefano Zampini   if (adapt) PetscValidHeaderSpecific(adapt, DM_CLASSID, 2);
403456cc5b7SMatthew G. Knepley   ierr = DMIsForest(dm, &isForest);CHKERRQ(ierr);
404456cc5b7SMatthew G. Knepley   if (!isForest) PetscFunctionReturn(0);
4051fd50544SStefano Zampini   if (adapt != NULL && dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
406ba936b91SToby Isaac   forest         = (DM_Forest*) dm->data;
407dffe73a3SToby Isaac   ierr           = DMForestGetAdaptivityForest(dm,&oldAdapt);CHKERRQ(ierr);
408193eb951SToby Isaac   adaptForest    = (DM_Forest*) (adapt ? adapt->data : NULL);
409193eb951SToby Isaac   oldAdaptForest = (DM_Forest*) (oldAdapt ? oldAdapt->data : NULL);
410dffe73a3SToby Isaac   if (adaptForest != oldAdaptForest) {
411dffe73a3SToby Isaac     ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr);
412dffe73a3SToby Isaac     ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr);
4131fd50544SStefano Zampini     if (forest->clearadaptivityforest) {ierr = (*forest->clearadaptivityforest)(dm);CHKERRQ(ierr);}
414dffe73a3SToby Isaac   }
41526d9498aSToby Isaac   switch (forest->adaptPurpose) {
416cd3c525cSToby Isaac   case DM_ADAPT_DETERMINE:
417ba936b91SToby Isaac     ierr          = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
418ba936b91SToby Isaac     ierr          = DMDestroy(&(forest->adapt));CHKERRQ(ierr);
419ba936b91SToby Isaac     forest->adapt = adapt;
42026d9498aSToby Isaac     break;
421a1b0c543SToby Isaac   case DM_ADAPT_REFINE:
42226d9498aSToby Isaac     ierr = DMSetCoarseDM(dm,adapt);CHKERRQ(ierr);
42326d9498aSToby Isaac     break;
424a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
425bf2d5fbbSStefano Zampini   case DM_ADAPT_COARSEN_LAST:
42626d9498aSToby Isaac     ierr = DMSetFineDM(dm,adapt);CHKERRQ(ierr);
42726d9498aSToby Isaac     break;
42826d9498aSToby Isaac   default:
42926d9498aSToby Isaac     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
43026d9498aSToby Isaac   }
431dd8e54a2SToby Isaac   PetscFunctionReturn(0);
432dd8e54a2SToby Isaac }
433dd8e54a2SToby Isaac 
4349be51f97SToby Isaac /*@
4359be51f97SToby Isaac   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.
4369be51f97SToby Isaac 
4379be51f97SToby Isaac   Not collective
4389be51f97SToby Isaac 
4399be51f97SToby Isaac   Input Parameter:
4409be51f97SToby Isaac . dm - the forest
4419be51f97SToby Isaac 
4429be51f97SToby Isaac   Output Parameter:
4439be51f97SToby Isaac . adapt - the forest from which dm is/was adapted
4449be51f97SToby Isaac 
4459be51f97SToby Isaac   Level: intermediate
4469be51f97SToby Isaac 
4479be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose()
4489be51f97SToby Isaac @*/
449ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
450dd8e54a2SToby Isaac {
451ba936b91SToby Isaac   DM_Forest      *forest;
45226d9498aSToby Isaac   PetscErrorCode ierr;
453dd8e54a2SToby Isaac 
454dd8e54a2SToby Isaac   PetscFunctionBegin;
455dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
456ba936b91SToby Isaac   forest = (DM_Forest*) dm->data;
45726d9498aSToby Isaac   switch (forest->adaptPurpose) {
458cd3c525cSToby Isaac   case DM_ADAPT_DETERMINE:
459ba936b91SToby Isaac     *adapt = forest->adapt;
46026d9498aSToby Isaac     break;
461a1b0c543SToby Isaac   case DM_ADAPT_REFINE:
46226d9498aSToby Isaac     ierr = DMGetCoarseDM(dm,adapt);CHKERRQ(ierr);
46326d9498aSToby Isaac     break;
464a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
465bf2d5fbbSStefano Zampini   case DM_ADAPT_COARSEN_LAST:
46626d9498aSToby Isaac     ierr = DMGetFineDM(dm,adapt);CHKERRQ(ierr);
46726d9498aSToby Isaac     break;
46826d9498aSToby Isaac   default:
46926d9498aSToby Isaac     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
47026d9498aSToby Isaac   }
47126d9498aSToby Isaac   PetscFunctionReturn(0);
47226d9498aSToby Isaac }
47326d9498aSToby Isaac 
4749be51f97SToby Isaac /*@
4759be51f97SToby Isaac   DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its
476a1b0c543SToby Isaac   source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening
477cd3c525cSToby Isaac   (DM_ADAPT_COARSEN), or undefined (DM_ADAPT_DETERMINE).  This only matters for the purposes of reference counting:
4789be51f97SToby Isaac   during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse
4799be51f97SToby Isaac   relationship (see DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does
4809be51f97SToby Isaac   not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can
4819be51f97SToby Isaac   cause a memory leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
4829be51f97SToby Isaac 
4839be51f97SToby Isaac   Logically collective on dm
4849be51f97SToby Isaac 
4859be51f97SToby Isaac   Input Parameters:
4869be51f97SToby Isaac + dm - the forest
487bf2d5fbbSStefano Zampini - purpose - the adaptivity purpose
4889be51f97SToby Isaac 
4899be51f97SToby Isaac   Level: advanced
4909be51f97SToby Isaac 
491bf2d5fbbSStefano Zampini .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest(), DMAdaptFlag
4929be51f97SToby Isaac @*/
493a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose)
49426d9498aSToby Isaac {
49526d9498aSToby Isaac   DM_Forest      *forest;
49626d9498aSToby Isaac   PetscErrorCode ierr;
49726d9498aSToby Isaac 
49826d9498aSToby Isaac   PetscFunctionBegin;
49926d9498aSToby Isaac   forest = (DM_Forest*) dm->data;
5009be51f97SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
50126d9498aSToby Isaac   if (purpose != forest->adaptPurpose) {
50226d9498aSToby Isaac     DM adapt;
50326d9498aSToby Isaac 
50426d9498aSToby Isaac     ierr = DMForestGetAdaptivityForest(dm,&adapt);CHKERRQ(ierr);
50526d9498aSToby Isaac     ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
50626d9498aSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
507f885a11aSToby Isaac 
50826d9498aSToby Isaac     forest->adaptPurpose = purpose;
509f885a11aSToby Isaac 
51026d9498aSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,adapt);CHKERRQ(ierr);
51126d9498aSToby Isaac     ierr = DMDestroy(&adapt);CHKERRQ(ierr);
51226d9498aSToby Isaac   }
513dd8e54a2SToby Isaac   PetscFunctionReturn(0);
514dd8e54a2SToby Isaac }
515dd8e54a2SToby Isaac 
51656c0450aSToby Isaac /*@
51756c0450aSToby Isaac   DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with
518bf2d5fbbSStefano Zampini   DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening (DM_ADAPT_COARSEN),
519bf2d5fbbSStefano Zampini   coarsening only the last level (DM_ADAPT_COARSEN_LAST) or undefined (DM_ADAPT_DETERMINE).
520bf2d5fbbSStefano Zampini   This only matters for the purposes of reference counting: during DMDestroy(), cyclic
52156c0450aSToby Isaac   references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see
52256c0450aSToby Isaac   DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does not maintain a
52356c0450aSToby Isaac   reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory
52456c0450aSToby Isaac   leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
52556c0450aSToby Isaac 
52656c0450aSToby Isaac   Not collective
52756c0450aSToby Isaac 
52856c0450aSToby Isaac   Input Parameter:
52956c0450aSToby Isaac . dm - the forest
53056c0450aSToby Isaac 
53156c0450aSToby Isaac   Output Parameter:
532bf2d5fbbSStefano Zampini . purpose - the adaptivity purpose
53356c0450aSToby Isaac 
53456c0450aSToby Isaac   Level: advanced
53556c0450aSToby Isaac 
536bf2d5fbbSStefano Zampini .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest(), DMAdaptFlag
53756c0450aSToby Isaac @*/
538a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose)
53956c0450aSToby Isaac {
54056c0450aSToby Isaac   DM_Forest *forest;
54156c0450aSToby Isaac 
54256c0450aSToby Isaac   PetscFunctionBegin;
54356c0450aSToby Isaac   forest   = (DM_Forest*) dm->data;
54456c0450aSToby Isaac   *purpose = forest->adaptPurpose;
54556c0450aSToby Isaac   PetscFunctionReturn(0);
54656c0450aSToby Isaac }
54756c0450aSToby Isaac 
5489be51f97SToby Isaac /*@
5499be51f97SToby Isaac   DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
5509be51f97SToby Isaac   cell adjacency (for the purposes of partitioning and overlap).
5519be51f97SToby Isaac 
5529be51f97SToby Isaac   Logically collective on dm
5539be51f97SToby Isaac 
5549be51f97SToby Isaac   Input Parameters:
5559be51f97SToby Isaac + dm - the forest
5569be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency)
5579be51f97SToby Isaac 
5589be51f97SToby Isaac   Level: intermediate
5599be51f97SToby Isaac 
5609be51f97SToby Isaac .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap()
5619be51f97SToby Isaac @*/
562dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
563dd8e54a2SToby Isaac {
564dd8e54a2SToby Isaac   PetscInt       dim;
565dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
566dd8e54a2SToby Isaac   PetscErrorCode ierr;
567dd8e54a2SToby Isaac 
568dd8e54a2SToby Isaac   PetscFunctionBegin;
569dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
570ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup");
571dd8e54a2SToby Isaac   if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim);
572dd8e54a2SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
573dd8e54a2SToby Isaac   if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim);
574dd8e54a2SToby Isaac   forest->adjDim = adjDim;
575dd8e54a2SToby Isaac   PetscFunctionReturn(0);
576dd8e54a2SToby Isaac }
577dd8e54a2SToby Isaac 
5789be51f97SToby Isaac /*@
5799be51f97SToby Isaac   DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that,
5809be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
5819be51f97SToby Isaac 
5829be51f97SToby Isaac   Logically collective on dm
5839be51f97SToby Isaac 
5849be51f97SToby Isaac   Input Parameters:
5859be51f97SToby Isaac + dm - the forest
5869be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
5879be51f97SToby Isaac 
5889be51f97SToby Isaac   Level: intermediate
5899be51f97SToby Isaac 
5909be51f97SToby Isaac .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension()
5919be51f97SToby Isaac @*/
592dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
593dd8e54a2SToby Isaac {
594dd8e54a2SToby Isaac   PetscInt       dim;
595dd8e54a2SToby Isaac   PetscErrorCode ierr;
596dd8e54a2SToby Isaac 
597dd8e54a2SToby Isaac   PetscFunctionBegin;
598dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
599dd8e54a2SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
600dd8e54a2SToby Isaac   ierr = DMForestSetAdjacencyDimension(dm,dim-adjCodim);CHKERRQ(ierr);
601dd8e54a2SToby Isaac   PetscFunctionReturn(0);
602dd8e54a2SToby Isaac }
603dd8e54a2SToby Isaac 
6049be51f97SToby Isaac /*@
6059be51f97SToby Isaac   DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
6069be51f97SToby Isaac   purposes of partitioning and overlap).
6079be51f97SToby Isaac 
6089be51f97SToby Isaac   Not collective
6099be51f97SToby Isaac 
6109be51f97SToby Isaac   Input Parameter:
6119be51f97SToby Isaac . dm - the forest
6129be51f97SToby Isaac 
6139be51f97SToby Isaac   Output Parameter:
6149be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency)
6159be51f97SToby Isaac 
6169be51f97SToby Isaac   Level: intermediate
6179be51f97SToby Isaac 
6189be51f97SToby Isaac .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap()
6199be51f97SToby Isaac @*/
620dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
621dd8e54a2SToby Isaac {
622dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
623dd8e54a2SToby Isaac 
624dd8e54a2SToby Isaac   PetscFunctionBegin;
625dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
626dd8e54a2SToby Isaac   PetscValidIntPointer(adjDim,2);
627dd8e54a2SToby Isaac   *adjDim = forest->adjDim;
628dd8e54a2SToby Isaac   PetscFunctionReturn(0);
629dd8e54a2SToby Isaac }
630dd8e54a2SToby Isaac 
6319be51f97SToby Isaac /*@
6329be51f97SToby Isaac   DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that,
6339be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
6349be51f97SToby Isaac 
6359be51f97SToby Isaac   Not collective
6369be51f97SToby Isaac 
6379be51f97SToby Isaac   Input Parameter:
6389be51f97SToby Isaac . dm - the forest
6399be51f97SToby Isaac 
6409be51f97SToby Isaac   Output Parameter:
6419be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
6429be51f97SToby Isaac 
6439be51f97SToby Isaac   Level: intermediate
6449be51f97SToby Isaac 
6459be51f97SToby Isaac .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension()
6469be51f97SToby Isaac @*/
647dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
648dd8e54a2SToby Isaac {
649dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
650dd8e54a2SToby Isaac   PetscInt       dim;
651dd8e54a2SToby Isaac   PetscErrorCode ierr;
652dd8e54a2SToby Isaac 
653dd8e54a2SToby Isaac   PetscFunctionBegin;
654dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
655dd8e54a2SToby Isaac   PetscValidIntPointer(adjCodim,2);
656dd8e54a2SToby Isaac   ierr      = DMGetDimension(dm,&dim);CHKERRQ(ierr);
657dd8e54a2SToby Isaac   *adjCodim = dim - forest->adjDim;
658dd8e54a2SToby Isaac   PetscFunctionReturn(0);
659dd8e54a2SToby Isaac }
660dd8e54a2SToby Isaac 
6619be51f97SToby Isaac /*@
6629be51f97SToby Isaac   DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel
6639be51f97SToby Isaac   partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding
6649be51f97SToby Isaac   adjacent cells
6659be51f97SToby Isaac 
6669be51f97SToby Isaac   Logically collective on dm
6679be51f97SToby Isaac 
6689be51f97SToby Isaac   Input Parameters:
6699be51f97SToby Isaac + dm - the forest
6709be51f97SToby Isaac - overlap - default 0
6719be51f97SToby Isaac 
6729be51f97SToby Isaac   Level: intermediate
6739be51f97SToby Isaac 
6749be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
6759be51f97SToby Isaac @*/
676dd8e54a2SToby Isaac PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap)
677dd8e54a2SToby Isaac {
678dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
679dd8e54a2SToby Isaac 
680dd8e54a2SToby Isaac   PetscFunctionBegin;
681dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
682ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup");
683dd8e54a2SToby Isaac   if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap);
684dd8e54a2SToby Isaac   forest->overlap = overlap;
685dd8e54a2SToby Isaac   PetscFunctionReturn(0);
686dd8e54a2SToby Isaac }
687dd8e54a2SToby Isaac 
6889be51f97SToby Isaac /*@
6899be51f97SToby Isaac   DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values
6909be51f97SToby Isaac   > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells
6919be51f97SToby Isaac 
6929be51f97SToby Isaac   Not collective
6939be51f97SToby Isaac 
6949be51f97SToby Isaac   Input Parameter:
6959be51f97SToby Isaac . dm - the forest
6969be51f97SToby Isaac 
6979be51f97SToby Isaac   Output Parameter:
6989be51f97SToby Isaac . overlap - default 0
6999be51f97SToby Isaac 
7009be51f97SToby Isaac   Level: intermediate
7019be51f97SToby Isaac 
7029be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
7039be51f97SToby Isaac @*/
704dd8e54a2SToby Isaac PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap)
705dd8e54a2SToby Isaac {
706dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
707dd8e54a2SToby Isaac 
708dd8e54a2SToby Isaac   PetscFunctionBegin;
709dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
710dd8e54a2SToby Isaac   PetscValidIntPointer(overlap,2);
711dd8e54a2SToby Isaac   *overlap = forest->overlap;
712dd8e54a2SToby Isaac   PetscFunctionReturn(0);
713dd8e54a2SToby Isaac }
714dd8e54a2SToby Isaac 
7159be51f97SToby Isaac /*@
7169be51f97SToby Isaac   DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base
7179be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest
7189be51f97SToby Isaac   (see DMForestGetAdaptivityForest()) this limits the amount of coarsening.
7199be51f97SToby Isaac 
7209be51f97SToby Isaac   Logically collective on dm
7219be51f97SToby Isaac 
7229be51f97SToby Isaac   Input Parameters:
7239be51f97SToby Isaac + dm - the forest
7249be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7259be51f97SToby Isaac 
7269be51f97SToby Isaac   Level: intermediate
7279be51f97SToby Isaac 
7289be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
7299be51f97SToby Isaac @*/
730dd8e54a2SToby Isaac PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement)
731dd8e54a2SToby Isaac {
732dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
733dd8e54a2SToby Isaac 
734dd8e54a2SToby Isaac   PetscFunctionBegin;
735dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
736ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup");
737dd8e54a2SToby Isaac   forest->minRefinement = minRefinement;
738dd8e54a2SToby Isaac   PetscFunctionReturn(0);
739dd8e54a2SToby Isaac }
740dd8e54a2SToby Isaac 
7419be51f97SToby Isaac /*@
7429be51f97SToby Isaac   DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see
7439be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest (see
7449be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of coarsening.
7459be51f97SToby Isaac 
7469be51f97SToby Isaac   Not collective
7479be51f97SToby Isaac 
7489be51f97SToby Isaac   Input Parameter:
7499be51f97SToby Isaac . dm - the forest
7509be51f97SToby Isaac 
7519be51f97SToby Isaac   Output Parameter:
7529be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7539be51f97SToby Isaac 
7549be51f97SToby Isaac   Level: intermediate
7559be51f97SToby Isaac 
7569be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
7579be51f97SToby Isaac @*/
758dd8e54a2SToby Isaac PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement)
759dd8e54a2SToby Isaac {
760dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
761dd8e54a2SToby Isaac 
762dd8e54a2SToby Isaac   PetscFunctionBegin;
763dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
764dd8e54a2SToby Isaac   PetscValidIntPointer(minRefinement,2);
765dd8e54a2SToby Isaac   *minRefinement = forest->minRefinement;
766dd8e54a2SToby Isaac   PetscFunctionReturn(0);
767dd8e54a2SToby Isaac }
768dd8e54a2SToby Isaac 
7699be51f97SToby Isaac /*@
7709be51f97SToby Isaac   DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base
7719be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.
7729be51f97SToby Isaac 
7739be51f97SToby Isaac   Logically collective on dm
7749be51f97SToby Isaac 
7759be51f97SToby Isaac   Input Parameters:
7769be51f97SToby Isaac + dm - the forest
7779be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7789be51f97SToby Isaac 
7799be51f97SToby Isaac   Level: intermediate
7809be51f97SToby Isaac 
7819be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
7829be51f97SToby Isaac @*/
78356ba9f64SToby Isaac PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement)
78456ba9f64SToby Isaac {
78556ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
78656ba9f64SToby Isaac 
78756ba9f64SToby Isaac   PetscFunctionBegin;
78856ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
78956ba9f64SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup");
79056ba9f64SToby Isaac   forest->initRefinement = initRefinement;
79156ba9f64SToby Isaac   PetscFunctionReturn(0);
79256ba9f64SToby Isaac }
79356ba9f64SToby Isaac 
7949be51f97SToby Isaac /*@
7959be51f97SToby Isaac   DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see
7969be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.
7979be51f97SToby Isaac 
7989be51f97SToby Isaac   Not collective
7999be51f97SToby Isaac 
8009be51f97SToby Isaac   Input Parameter:
8019be51f97SToby Isaac . dm - the forest
8029be51f97SToby Isaac 
80301d2d390SJose E. Roman   Output Parameter:
804bf2d5fbbSStefano Zampini . initRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8059be51f97SToby Isaac 
8069be51f97SToby Isaac   Level: intermediate
8079be51f97SToby Isaac 
8089be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
8099be51f97SToby Isaac @*/
81056ba9f64SToby Isaac PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement)
81156ba9f64SToby Isaac {
81256ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
81356ba9f64SToby Isaac 
81456ba9f64SToby Isaac   PetscFunctionBegin;
81556ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
81656ba9f64SToby Isaac   PetscValidIntPointer(initRefinement,2);
81756ba9f64SToby Isaac   *initRefinement = forest->initRefinement;
81856ba9f64SToby Isaac   PetscFunctionReturn(0);
81956ba9f64SToby Isaac }
82056ba9f64SToby Isaac 
8219be51f97SToby Isaac /*@
8229be51f97SToby Isaac   DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base
8239be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest
8249be51f97SToby Isaac   (see DMForestGetAdaptivityForest()), this limits the amount of refinement.
8259be51f97SToby Isaac 
8269be51f97SToby Isaac   Logically collective on dm
8279be51f97SToby Isaac 
8289be51f97SToby Isaac   Input Parameters:
8299be51f97SToby Isaac + dm - the forest
8309be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8319be51f97SToby Isaac 
8329be51f97SToby Isaac   Level: intermediate
8339be51f97SToby Isaac 
8349be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM()
8359be51f97SToby Isaac @*/
836c7eeac06SToby Isaac PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement)
837dd8e54a2SToby Isaac {
838dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
839dd8e54a2SToby Isaac 
840dd8e54a2SToby Isaac   PetscFunctionBegin;
841dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
842ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup");
843c7eeac06SToby Isaac   forest->maxRefinement = maxRefinement;
844dd8e54a2SToby Isaac   PetscFunctionReturn(0);
845dd8e54a2SToby Isaac }
846dd8e54a2SToby Isaac 
8479be51f97SToby Isaac /*@
8489be51f97SToby Isaac   DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see
8499be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest (see
8509be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of refinement.
8519be51f97SToby Isaac 
8529be51f97SToby Isaac   Not collective
8539be51f97SToby Isaac 
8549be51f97SToby Isaac   Input Parameter:
8559be51f97SToby Isaac . dm - the forest
8569be51f97SToby Isaac 
8579be51f97SToby Isaac   Output Parameter:
8589be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8599be51f97SToby Isaac 
8609be51f97SToby Isaac   Level: intermediate
8619be51f97SToby Isaac 
8629be51f97SToby Isaac .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
8639be51f97SToby Isaac @*/
864c7eeac06SToby Isaac PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement)
865dd8e54a2SToby Isaac {
866dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
867dd8e54a2SToby Isaac 
868dd8e54a2SToby Isaac   PetscFunctionBegin;
869dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
870c7eeac06SToby Isaac   PetscValidIntPointer(maxRefinement,2);
871c7eeac06SToby Isaac   *maxRefinement = forest->maxRefinement;
872dd8e54a2SToby Isaac   PetscFunctionReturn(0);
873dd8e54a2SToby Isaac }
874c7eeac06SToby Isaac 
8759be51f97SToby Isaac /*@C
8769be51f97SToby Isaac   DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes.
8779be51f97SToby Isaac   Subtypes of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree
8789be51f97SToby Isaac   for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to
8799be51f97SToby Isaac   specify refinement/coarsening.
8809be51f97SToby Isaac 
8819be51f97SToby Isaac   Logically collective on dm
8829be51f97SToby Isaac 
8839be51f97SToby Isaac   Input Parameters:
8849be51f97SToby Isaac + dm - the forest
8859be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL
8869be51f97SToby Isaac 
8879be51f97SToby Isaac   Level: advanced
8889be51f97SToby Isaac 
8899be51f97SToby Isaac .seealso: DMForestGetAdaptivityStrategy()
8909be51f97SToby Isaac @*/
891c7eeac06SToby Isaac PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy)
892c7eeac06SToby Isaac {
893c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
894c7eeac06SToby Isaac   PetscErrorCode ierr;
895c7eeac06SToby Isaac 
896c7eeac06SToby Isaac   PetscFunctionBegin;
897c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
898c7eeac06SToby Isaac   ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr);
899a73e2921SToby Isaac   ierr = PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy);CHKERRQ(ierr);
900c7eeac06SToby Isaac   PetscFunctionReturn(0);
901c7eeac06SToby Isaac }
902c7eeac06SToby Isaac 
9039be51f97SToby Isaac /*@C
9049be51f97SToby Isaac   DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes.  Subtypes
9059be51f97SToby Isaac   of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all
9069be51f97SToby Isaac   processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only
9079be51f97SToby Isaac   one process needs to specify refinement/coarsening.
9089be51f97SToby Isaac 
9099be51f97SToby Isaac   Not collective
9109be51f97SToby Isaac 
9119be51f97SToby Isaac   Input Parameter:
9129be51f97SToby Isaac . dm - the forest
9139be51f97SToby Isaac 
9149be51f97SToby Isaac   Output Parameter:
9159be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL)
9169be51f97SToby Isaac 
9179be51f97SToby Isaac   Level: advanced
9189be51f97SToby Isaac 
9199be51f97SToby Isaac .seealso: DMForestSetAdaptivityStrategy()
9209be51f97SToby Isaac @*/
921c7eeac06SToby Isaac PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
922c7eeac06SToby Isaac {
923c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
924c7eeac06SToby Isaac 
925c7eeac06SToby Isaac   PetscFunctionBegin;
926c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
927c7eeac06SToby Isaac   PetscValidPointer(adaptStrategy,2);
928c7eeac06SToby Isaac   *adaptStrategy = forest->adaptStrategy;
929c7eeac06SToby Isaac   PetscFunctionReturn(0);
930c7eeac06SToby Isaac }
931c7eeac06SToby Isaac 
9322a133e43SToby Isaac /*@
9332a133e43SToby Isaac   DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning,
9342a133e43SToby Isaac   etc.) was successful.  PETSC_FALSE indicates that the post-adaptation forest is the same as the pre-adpatation
9352a133e43SToby Isaac   forest.  A requested adaptation may have been unsuccessful if, for example, the requested refinement would have
9362a133e43SToby Isaac   exceeded the maximum refinement level.
9372a133e43SToby Isaac 
9382a133e43SToby Isaac   Collective on dm
9392a133e43SToby Isaac 
9402a133e43SToby Isaac   Input Parameter:
9412a133e43SToby Isaac 
9422a133e43SToby Isaac . dm - the post-adaptation forest
9432a133e43SToby Isaac 
9442a133e43SToby Isaac   Output Parameter:
9452a133e43SToby Isaac 
9462a133e43SToby Isaac . success - PETSC_TRUE if the post-adaptation forest is different from the pre-adaptation forest.
9472a133e43SToby Isaac 
9482a133e43SToby Isaac   Level: intermediate
9492a133e43SToby Isaac 
9502a133e43SToby Isaac .see
9512a133e43SToby Isaac @*/
9522a133e43SToby Isaac PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success)
9532a133e43SToby Isaac {
9542a133e43SToby Isaac   DM_Forest      *forest;
9552a133e43SToby Isaac   PetscErrorCode ierr;
9562a133e43SToby Isaac 
9572a133e43SToby Isaac   PetscFunctionBegin;
9582a133e43SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9592a133e43SToby Isaac   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() has not been called yet.");
9602a133e43SToby Isaac   forest = (DM_Forest *) dm->data;
9612a133e43SToby Isaac   ierr = (forest->getadaptivitysuccess)(dm,success);CHKERRQ(ierr);
9622a133e43SToby Isaac   PetscFunctionReturn(0);
9632a133e43SToby Isaac }
9642a133e43SToby Isaac 
965bf9b5d84SToby Isaac /*@
966bf9b5d84SToby Isaac   DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed
967bf9b5d84SToby 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().
968bf9b5d84SToby Isaac 
969bf9b5d84SToby Isaac   Logically collective on dm
970bf9b5d84SToby Isaac 
971bf9b5d84SToby Isaac   Input Parameters:
972bf9b5d84SToby Isaac + dm - the post-adaptation forest
973bf9b5d84SToby Isaac - computeSF - default PETSC_TRUE
974bf9b5d84SToby Isaac 
975bf9b5d84SToby Isaac   Level: advanced
976bf9b5d84SToby Isaac 
977bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
978bf9b5d84SToby Isaac @*/
979bf9b5d84SToby Isaac PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
980bf9b5d84SToby Isaac {
981bf9b5d84SToby Isaac   DM_Forest *forest;
982bf9b5d84SToby Isaac 
983bf9b5d84SToby Isaac   PetscFunctionBegin;
984bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
985bf9b5d84SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called");
986bf9b5d84SToby Isaac   forest                 = (DM_Forest*) dm->data;
987bf9b5d84SToby Isaac   forest->computeAdaptSF = computeSF;
988bf9b5d84SToby Isaac   PetscFunctionReturn(0);
989bf9b5d84SToby Isaac }
990bf9b5d84SToby Isaac 
9910eb7e1eaSToby Isaac PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time)
99280b27e07SToby Isaac {
99380b27e07SToby Isaac   DM_Forest      *forest;
99480b27e07SToby Isaac   PetscErrorCode ierr;
99580b27e07SToby Isaac 
99680b27e07SToby Isaac   PetscFunctionBegin;
99780b27e07SToby Isaac   PetscValidHeaderSpecific(dmIn   ,DM_CLASSID  ,1);
99880b27e07SToby Isaac   PetscValidHeaderSpecific(vecIn  ,VEC_CLASSID ,2);
99980b27e07SToby Isaac   PetscValidHeaderSpecific(dmOut  ,DM_CLASSID  ,3);
100080b27e07SToby Isaac   PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,4);
100180b27e07SToby Isaac   forest = (DM_Forest *) dmIn->data;
100280b27e07SToby Isaac   if (!forest->transfervec) SETERRQ(PetscObjectComm((PetscObject)dmIn),PETSC_ERR_SUP,"DMForestTransferVec() not implemented");
10030eb7e1eaSToby Isaac   ierr = (forest->transfervec)(dmIn,vecIn,dmOut,vecOut,useBCs,time);CHKERRQ(ierr);
100480b27e07SToby Isaac   PetscFunctionReturn(0);
100580b27e07SToby Isaac }
100680b27e07SToby Isaac 
1007ac34a06fSStefano Zampini PetscErrorCode DMForestTransferVecFromBase(DM dm, Vec vecIn, Vec vecOut)
1008ac34a06fSStefano Zampini {
1009ac34a06fSStefano Zampini   DM_Forest      *forest;
1010ac34a06fSStefano Zampini   PetscErrorCode ierr;
1011ac34a06fSStefano Zampini 
1012ac34a06fSStefano Zampini   PetscFunctionBegin;
1013ac34a06fSStefano Zampini   PetscValidHeaderSpecific(dm   ,DM_CLASSID  ,1);
1014ac34a06fSStefano Zampini   PetscValidHeaderSpecific(vecIn  ,VEC_CLASSID ,2);
1015ac34a06fSStefano Zampini   PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,3);
1016ac34a06fSStefano Zampini   forest = (DM_Forest *) dm->data;
1017ac34a06fSStefano Zampini   if (!forest->transfervecfrombase) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DMForestTransferVecFromBase() not implemented");
1018ac34a06fSStefano Zampini   ierr = (forest->transfervecfrombase)(dm,vecIn,vecOut);CHKERRQ(ierr);
1019ac34a06fSStefano Zampini   PetscFunctionReturn(0);
1020ac34a06fSStefano Zampini }
1021ac34a06fSStefano Zampini 
1022bf9b5d84SToby Isaac /*@
1023bf9b5d84SToby Isaac   DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the
1024bf9b5d84SToby Isaac   pre-adaptation forest to the post-adaptiation forest.  After DMSetUp() is called, these transfer PetscSFs can be
1025bf9b5d84SToby Isaac   accessed with DMForestGetAdaptivitySF().
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:
1033bf9b5d84SToby Isaac . computeSF - default PETSC_TRUE
1034bf9b5d84SToby Isaac 
1035bf9b5d84SToby Isaac   Level: advanced
1036bf9b5d84SToby Isaac 
1037bf9b5d84SToby Isaac .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
1038bf9b5d84SToby Isaac @*/
1039bf9b5d84SToby Isaac PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1040bf9b5d84SToby Isaac {
1041bf9b5d84SToby Isaac   DM_Forest *forest;
1042bf9b5d84SToby Isaac 
1043bf9b5d84SToby Isaac   PetscFunctionBegin;
1044bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1045bf9b5d84SToby Isaac   forest     = (DM_Forest*) dm->data;
1046bf9b5d84SToby Isaac   *computeSF = forest->computeAdaptSF;
1047bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1048bf9b5d84SToby Isaac }
1049bf9b5d84SToby Isaac 
1050bf9b5d84SToby Isaac /*@
1051bf9b5d84SToby Isaac   DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest.
1052bf9b5d84SToby Isaac   Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be
1053bf9b5d84SToby Isaac   some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa.  Therefore there are two
1054bf9b5d84SToby Isaac   PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates
1055bf9b5d84SToby Isaac   pre-adaptation fine cells to post-adaptation coarse cells.
1056bf9b5d84SToby Isaac 
1057bf9b5d84SToby Isaac   Not collective
1058bf9b5d84SToby Isaac 
1059bf9b5d84SToby Isaac   Input Parameter:
1060bf9b5d84SToby Isaac   dm - the post-adaptation forest
1061bf9b5d84SToby Isaac 
1062bf9b5d84SToby Isaac   Output Parameter:
10630f17b9e3SToby Isaac   preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
10640f17b9e3SToby Isaac   coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-
1065bf9b5d84SToby Isaac 
1066bf9b5d84SToby Isaac   Level: advanced
1067bf9b5d84SToby Isaac 
1068bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF()
1069bf9b5d84SToby Isaac @*/
10700f17b9e3SToby Isaac PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1071bf9b5d84SToby Isaac {
1072bf9b5d84SToby Isaac   DM_Forest      *forest;
1073bf9b5d84SToby Isaac   PetscErrorCode ierr;
1074bf9b5d84SToby Isaac 
1075bf9b5d84SToby Isaac   PetscFunctionBegin;
1076bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1077bf9b5d84SToby Isaac   ierr   = DMSetUp(dm);CHKERRQ(ierr);
1078bf9b5d84SToby Isaac   forest = (DM_Forest*) dm->data;
1079f885a11aSToby Isaac   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1080f885a11aSToby Isaac   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
1081bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1082bf9b5d84SToby Isaac }
1083bf9b5d84SToby Isaac 
10849be51f97SToby Isaac /*@
10859be51f97SToby Isaac   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
10869be51f97SToby Isaac   indicate that the diameter of neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may
10879be51f97SToby Isaac   only support one particular choice of grading factor.
10889be51f97SToby Isaac 
10899be51f97SToby Isaac   Logically collective on dm
10909be51f97SToby Isaac 
10919be51f97SToby Isaac   Input Parameters:
10929be51f97SToby Isaac + dm - the forest
10939be51f97SToby Isaac - grade - the grading factor
10949be51f97SToby Isaac 
10959be51f97SToby Isaac   Level: advanced
10969be51f97SToby Isaac 
10979be51f97SToby Isaac .seealso: DMForestGetGradeFactor()
10989be51f97SToby Isaac @*/
1099c7eeac06SToby Isaac PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1100c7eeac06SToby Isaac {
1101c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1102c7eeac06SToby Isaac 
1103c7eeac06SToby Isaac   PetscFunctionBegin;
1104c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1105ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup");
1106c7eeac06SToby Isaac   forest->gradeFactor = grade;
1107c7eeac06SToby Isaac   PetscFunctionReturn(0);
1108c7eeac06SToby Isaac }
1109c7eeac06SToby Isaac 
11109be51f97SToby Isaac /*@
11119be51f97SToby Isaac   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
11129be51f97SToby Isaac   neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may only support one particular
11139be51f97SToby Isaac   choice of grading factor.
11149be51f97SToby Isaac 
11159be51f97SToby Isaac   Not collective
11169be51f97SToby Isaac 
11179be51f97SToby Isaac   Input Parameter:
11189be51f97SToby Isaac . dm - the forest
11199be51f97SToby Isaac 
11209be51f97SToby Isaac   Output Parameter:
11219be51f97SToby Isaac . grade - the grading factor
11229be51f97SToby Isaac 
11239be51f97SToby Isaac   Level: advanced
11249be51f97SToby Isaac 
11259be51f97SToby Isaac .seealso: DMForestSetGradeFactor()
11269be51f97SToby Isaac @*/
1127c7eeac06SToby Isaac PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1128c7eeac06SToby Isaac {
1129c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1130c7eeac06SToby Isaac 
1131c7eeac06SToby Isaac   PetscFunctionBegin;
1132c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1133c7eeac06SToby Isaac   PetscValidIntPointer(grade,2);
1134c7eeac06SToby Isaac   *grade = forest->gradeFactor;
1135c7eeac06SToby Isaac   PetscFunctionReturn(0);
1136c7eeac06SToby Isaac }
1137c7eeac06SToby Isaac 
11389be51f97SToby Isaac /*@
11399be51f97SToby Isaac   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
11409be51f97SToby Isaac   the cell weight (see DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be
11419be51f97SToby Isaac   (cellWeight) * (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on
11429be51f97SToby Isaac   its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
11439be51f97SToby Isaac   computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11449be51f97SToby Isaac 
11459be51f97SToby Isaac   Logically collective on dm
11469be51f97SToby Isaac 
11479be51f97SToby Isaac   Input Parameters:
11489be51f97SToby Isaac + dm - the forest
11499be51f97SToby Isaac - weightsFactors - default 1.
11509be51f97SToby Isaac 
11519be51f97SToby Isaac   Level: advanced
11529be51f97SToby Isaac 
11539be51f97SToby Isaac .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights()
11549be51f97SToby Isaac @*/
1155ef51cf95SToby Isaac PetscErrorCode DMForestSetCellWeightFactor(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);
1161ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup");
1162c7eeac06SToby Isaac   forest->weightsFactor = weightsFactor;
1163c7eeac06SToby Isaac   PetscFunctionReturn(0);
1164c7eeac06SToby Isaac }
1165c7eeac06SToby Isaac 
11669be51f97SToby Isaac /*@
11679be51f97SToby Isaac   DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
11689be51f97SToby Isaac   DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be (cellWeight) *
11699be51f97SToby Isaac   (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on its level; a
11709be51f97SToby Isaac   factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation
11719be51f97SToby Isaac   associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11729be51f97SToby Isaac 
11739be51f97SToby Isaac   Not collective
11749be51f97SToby Isaac 
11759be51f97SToby Isaac   Input Parameter:
11769be51f97SToby Isaac . dm - the forest
11779be51f97SToby Isaac 
11789be51f97SToby Isaac   Output Parameter:
11799be51f97SToby Isaac . weightsFactors - default 1.
11809be51f97SToby Isaac 
11819be51f97SToby Isaac   Level: advanced
11829be51f97SToby Isaac 
11839be51f97SToby Isaac .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights()
11849be51f97SToby Isaac @*/
1185ef51cf95SToby Isaac PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1186c7eeac06SToby Isaac {
1187c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1188c7eeac06SToby Isaac 
1189c7eeac06SToby Isaac   PetscFunctionBegin;
1190c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1191c7eeac06SToby Isaac   PetscValidRealPointer(weightsFactor,2);
1192c7eeac06SToby Isaac   *weightsFactor = forest->weightsFactor;
1193c7eeac06SToby Isaac   PetscFunctionReturn(0);
1194c7eeac06SToby Isaac }
1195c7eeac06SToby Isaac 
11969be51f97SToby Isaac /*@
11979be51f97SToby Isaac   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
11989be51f97SToby Isaac 
11999be51f97SToby Isaac   Not collective
12009be51f97SToby Isaac 
12019be51f97SToby Isaac   Input Parameter:
12029be51f97SToby Isaac . dm - the forest
12039be51f97SToby Isaac 
12049be51f97SToby Isaac   Output Parameters:
12059be51f97SToby Isaac + cStart - the first cell on this process
12069be51f97SToby Isaac - cEnd - one after the final cell on this process
12079be51f97SToby Isaac 
12081a244344SSatish Balay   Level: intermediate
12099be51f97SToby Isaac 
12109be51f97SToby Isaac .seealso: DMForestGetCellSF()
12119be51f97SToby Isaac @*/
1212c7eeac06SToby Isaac PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1213c7eeac06SToby Isaac {
1214c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1215c7eeac06SToby Isaac   PetscErrorCode ierr;
1216c7eeac06SToby Isaac 
1217c7eeac06SToby Isaac   PetscFunctionBegin;
1218c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1219c7eeac06SToby Isaac   PetscValidIntPointer(cStart,2);
1220064a246eSJacob Faibussowitsch   PetscValidIntPointer(cEnd,3);
1221c7eeac06SToby Isaac   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) {
1222c7eeac06SToby Isaac     ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr);
1223c7eeac06SToby Isaac   }
1224c7eeac06SToby Isaac   *cStart =  forest->cStart;
1225c7eeac06SToby Isaac   *cEnd   =  forest->cEnd;
1226c7eeac06SToby Isaac   PetscFunctionReturn(0);
1227c7eeac06SToby Isaac }
1228c7eeac06SToby Isaac 
12299be51f97SToby Isaac /*@
12309be51f97SToby Isaac   DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes
12319be51f97SToby Isaac 
12329be51f97SToby Isaac   Not collective
12339be51f97SToby Isaac 
12349be51f97SToby Isaac   Input Parameter:
12359be51f97SToby Isaac . dm - the forest
12369be51f97SToby Isaac 
12379be51f97SToby Isaac   Output Parameter:
12389be51f97SToby Isaac . cellSF - the PetscSF
12399be51f97SToby Isaac 
12401a244344SSatish Balay   Level: intermediate
12419be51f97SToby Isaac 
12429be51f97SToby Isaac .seealso: DMForestGetCellChart()
12439be51f97SToby Isaac @*/
1244c7eeac06SToby Isaac PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1245c7eeac06SToby Isaac {
1246c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1247c7eeac06SToby Isaac   PetscErrorCode ierr;
1248c7eeac06SToby Isaac 
1249c7eeac06SToby Isaac   PetscFunctionBegin;
1250c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1251c7eeac06SToby Isaac   PetscValidPointer(cellSF,2);
1252c7eeac06SToby Isaac   if ((!forest->cellSF) && forest->createcellsf) {
1253c7eeac06SToby Isaac     ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr);
1254c7eeac06SToby Isaac   }
1255c7eeac06SToby Isaac   *cellSF = forest->cellSF;
1256c7eeac06SToby Isaac   PetscFunctionReturn(0);
1257c7eeac06SToby Isaac }
1258c7eeac06SToby Isaac 
12599be51f97SToby Isaac /*@C
12609be51f97SToby Isaac   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
12619be51f97SToby Isaac   DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination).  The
1262cd3c525cSToby Isaac   interpretation of the label values is up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP,
1263cd3c525cSToby Isaac   DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have been reserved as choices that should be accepted by all subtypes.
12649be51f97SToby Isaac 
12659be51f97SToby Isaac   Logically collective on dm
12669be51f97SToby Isaac 
12679be51f97SToby Isaac   Input Parameters:
12689be51f97SToby Isaac - dm - the forest
1269a1b0c543SToby Isaac + adaptLabel - the label in the pre-adaptation forest
12709be51f97SToby Isaac 
12719be51f97SToby Isaac   Level: intermediate
12729be51f97SToby Isaac 
12739be51f97SToby Isaac .seealso DMForestGetAdaptivityLabel()
12749be51f97SToby Isaac @*/
1275a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel)
1276c7eeac06SToby Isaac {
1277c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1278c7eeac06SToby Isaac   PetscErrorCode ierr;
1279c7eeac06SToby Isaac 
1280c7eeac06SToby Isaac   PetscFunctionBegin;
1281c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
12821fd50544SStefano Zampini   if (adaptLabel) PetscValidHeaderSpecific(adaptLabel,DMLABEL_CLASSID,2);
1283d67d17b1SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject)adaptLabel);CHKERRQ(ierr);
12841fd50544SStefano Zampini   ierr = DMLabelDestroy(&forest->adaptLabel);CHKERRQ(ierr);
12851fd50544SStefano Zampini   forest->adaptLabel = adaptLabel;
1286c7eeac06SToby Isaac   PetscFunctionReturn(0);
1287c7eeac06SToby Isaac }
1288c7eeac06SToby Isaac 
12899be51f97SToby Isaac /*@C
12909be51f97SToby Isaac   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that
12919be51f97SToby Isaac   holds the adaptation flags (refinement, coarsening, or some combination).  The interpretation of the label values is
1292cd3c525cSToby Isaac   up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have
1293cd3c525cSToby Isaac   been reserved as choices that should be accepted by all subtypes.
12949be51f97SToby Isaac 
12959be51f97SToby Isaac   Not collective
12969be51f97SToby Isaac 
12979be51f97SToby Isaac   Input Parameter:
12989be51f97SToby Isaac . dm - the forest
12999be51f97SToby Isaac 
13009be51f97SToby Isaac   Output Parameter:
13019be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest
13029be51f97SToby Isaac 
13039be51f97SToby Isaac   Level: intermediate
13049be51f97SToby Isaac 
13059be51f97SToby Isaac .seealso DMForestSetAdaptivityLabel()
13069be51f97SToby Isaac @*/
1307a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel)
1308c7eeac06SToby Isaac {
1309c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1310c7eeac06SToby Isaac 
1311c7eeac06SToby Isaac   PetscFunctionBegin;
1312c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1313ba936b91SToby Isaac   *adaptLabel = forest->adaptLabel;
1314c7eeac06SToby Isaac   PetscFunctionReturn(0);
1315c7eeac06SToby Isaac }
1316c7eeac06SToby Isaac 
13179be51f97SToby Isaac /*@
13189be51f97SToby Isaac   DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
13199be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
13209be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
13219be51f97SToby Isaac   of 1.
13229be51f97SToby Isaac 
13239be51f97SToby Isaac   Logically collective on dm
13249be51f97SToby Isaac 
13259be51f97SToby Isaac   Input Parameters:
13269be51f97SToby Isaac + dm - the forest
13279be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
13289be51f97SToby Isaac - copyMode - how weights should reference weights
13299be51f97SToby Isaac 
13309be51f97SToby Isaac   Level: advanced
13319be51f97SToby Isaac 
13329be51f97SToby Isaac .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity()
13339be51f97SToby Isaac @*/
1334c7eeac06SToby Isaac PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1335c7eeac06SToby Isaac {
1336c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1337c7eeac06SToby Isaac   PetscInt       cStart, cEnd;
1338c7eeac06SToby Isaac   PetscErrorCode ierr;
1339c7eeac06SToby Isaac 
1340c7eeac06SToby Isaac   PetscFunctionBegin;
1341c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1342c7eeac06SToby Isaac   ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr);
1343c7eeac06SToby Isaac   if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd);
1344c7eeac06SToby Isaac   if (copyMode == PETSC_COPY_VALUES) {
1345c7eeac06SToby Isaac     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) {
1346c7eeac06SToby Isaac       ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr);
1347c7eeac06SToby Isaac     }
1348580bdb30SBarry Smith     ierr                        = PetscArraycpy(forest->cellWeights,weights,cEnd-cStart);CHKERRQ(ierr);
1349c7eeac06SToby Isaac     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
1350c7eeac06SToby Isaac     PetscFunctionReturn(0);
1351c7eeac06SToby Isaac   }
1352c7eeac06SToby Isaac   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) {
1353c7eeac06SToby Isaac     ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr);
1354c7eeac06SToby Isaac   }
1355c7eeac06SToby Isaac   forest->cellWeights         = weights;
1356c7eeac06SToby Isaac   forest->cellWeightsCopyMode = copyMode;
1357c7eeac06SToby Isaac   PetscFunctionReturn(0);
1358c7eeac06SToby Isaac }
1359c7eeac06SToby Isaac 
13609be51f97SToby Isaac /*@
13619be51f97SToby Isaac   DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
13629be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
13639be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
13649be51f97SToby Isaac   of 1.
13659be51f97SToby Isaac 
13669be51f97SToby Isaac   Not collective
13679be51f97SToby Isaac 
13689be51f97SToby Isaac   Input Parameter:
13699be51f97SToby Isaac . dm - the forest
13709be51f97SToby Isaac 
13719be51f97SToby Isaac   Output Parameter:
13729be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
13739be51f97SToby Isaac 
13749be51f97SToby Isaac   Level: advanced
13759be51f97SToby Isaac 
13769be51f97SToby Isaac .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity()
13779be51f97SToby Isaac @*/
1378c7eeac06SToby Isaac PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1379c7eeac06SToby Isaac {
1380c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1381c7eeac06SToby Isaac 
1382c7eeac06SToby Isaac   PetscFunctionBegin;
1383c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1384c7eeac06SToby Isaac   PetscValidPointer(weights,2);
1385c7eeac06SToby Isaac   *weights = forest->cellWeights;
1386c7eeac06SToby Isaac   PetscFunctionReturn(0);
1387c7eeac06SToby Isaac }
1388c7eeac06SToby Isaac 
13899be51f97SToby Isaac /*@
13909be51f97SToby Isaac   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
13919be51f97SToby Isaac   a pre-adaptation forest (see DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each
13929be51f97SToby Isaac   process's cells to the process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that
13939be51f97SToby Isaac   the current process should not have any cells after repartitioning.
13949be51f97SToby Isaac 
13959be51f97SToby Isaac   Logically Collective on dm
13969be51f97SToby Isaac 
13979be51f97SToby Isaac   Input parameters:
13989be51f97SToby Isaac + dm - the forest
13999be51f97SToby Isaac - capacity - this process's capacity
14009be51f97SToby Isaac 
14019be51f97SToby Isaac   Level: advanced
14029be51f97SToby Isaac 
14039be51f97SToby Isaac .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
14049be51f97SToby Isaac @*/
1405c7eeac06SToby Isaac PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1406c7eeac06SToby Isaac {
1407c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1408c7eeac06SToby Isaac 
1409c7eeac06SToby Isaac   PetscFunctionBegin;
1410c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1411ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup");
1412c7eeac06SToby Isaac   if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity);
1413c7eeac06SToby Isaac   forest->weightCapacity = capacity;
1414c7eeac06SToby Isaac   PetscFunctionReturn(0);
1415c7eeac06SToby Isaac }
1416c7eeac06SToby Isaac 
14179be51f97SToby Isaac /*@
14189be51f97SToby Isaac   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
14199be51f97SToby Isaac   DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each process's cells to the
14209be51f97SToby Isaac   process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that the current process
14219be51f97SToby Isaac   should not have any cells after repartitioning.
14229be51f97SToby Isaac 
14239be51f97SToby Isaac   Not collective
14249be51f97SToby Isaac 
14259be51f97SToby Isaac   Input parameter:
14269be51f97SToby Isaac . dm - the forest
14279be51f97SToby Isaac 
14289be51f97SToby Isaac   Output parameter:
14299be51f97SToby Isaac . capacity - this process's capacity
14309be51f97SToby Isaac 
14319be51f97SToby Isaac   Level: advanced
14329be51f97SToby Isaac 
14339be51f97SToby Isaac .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
14349be51f97SToby Isaac @*/
1435c7eeac06SToby Isaac PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1436c7eeac06SToby Isaac {
1437c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1438c7eeac06SToby Isaac 
1439c7eeac06SToby Isaac   PetscFunctionBegin;
1440c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1441c7eeac06SToby Isaac   PetscValidRealPointer(capacity,2);
1442c7eeac06SToby Isaac   *capacity = forest->weightCapacity;
1443c7eeac06SToby Isaac   PetscFunctionReturn(0);
1444c7eeac06SToby Isaac }
1445c7eeac06SToby Isaac 
14460709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm)
1447db4d5e8cSToby Isaac {
144856ba9f64SToby Isaac   PetscBool                  flg, flg1, flg2, flg3, flg4;
1449dd8e54a2SToby Isaac   DMForestTopology           oldTopo;
1450c7eeac06SToby Isaac   char                       stringBuffer[256];
1451dd8e54a2SToby Isaac   PetscViewer                viewer;
1452dd8e54a2SToby Isaac   PetscViewerFormat          format;
145356ba9f64SToby Isaac   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1454c7eeac06SToby Isaac   PetscReal                  weightsFactor;
1455c7eeac06SToby Isaac   DMForestAdaptivityStrategy adaptStrategy;
1456db4d5e8cSToby Isaac   PetscErrorCode             ierr;
1457db4d5e8cSToby Isaac 
1458db4d5e8cSToby Isaac   PetscFunctionBegin;
1459db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1460dd8e54a2SToby Isaac   ierr                         = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr);
1461a3eda1eaSToby Isaac   ierr                         = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr);
1462589a23caSBarry Smith   ierr                         = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,sizeof(stringBuffer),&flg1);CHKERRQ(ierr);
146356ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr);
146456ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr);
146556ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr);
1466f885a11aSToby 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}");
146756ba9f64SToby Isaac   if (flg1) {
146856ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr);
146956ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
147020e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
147156ba9f64SToby Isaac   }
147256ba9f64SToby Isaac   if (flg2) {
1473dd8e54a2SToby Isaac     DM base;
1474dd8e54a2SToby Isaac 
1475dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr);
1476dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1477dd8e54a2SToby Isaac     ierr = DMLoad(base,viewer);CHKERRQ(ierr);
1478dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1479dd8e54a2SToby Isaac     ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr);
1480dd8e54a2SToby Isaac     ierr = DMDestroy(&base);CHKERRQ(ierr);
148156ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
148220e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
1483dd8e54a2SToby Isaac   }
148456ba9f64SToby Isaac   if (flg3) {
1485dd8e54a2SToby Isaac     DM coarse;
1486dd8e54a2SToby Isaac 
1487dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr);
1488dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1489dd8e54a2SToby Isaac     ierr = DMLoad(coarse,viewer);CHKERRQ(ierr);
1490dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
149120e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr);
1492dd8e54a2SToby Isaac     ierr = DMDestroy(&coarse);CHKERRQ(ierr);
149356ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
149456ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1495dd8e54a2SToby Isaac   }
149656ba9f64SToby Isaac   if (flg4) {
1497dd8e54a2SToby Isaac     DM fine;
1498dd8e54a2SToby Isaac 
1499dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr);
1500dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1501dd8e54a2SToby Isaac     ierr = DMLoad(fine,viewer);CHKERRQ(ierr);
1502dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
150320e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr);
1504dd8e54a2SToby Isaac     ierr = DMDestroy(&fine);CHKERRQ(ierr);
150556ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
150656ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1507dd8e54a2SToby Isaac   }
1508dd8e54a2SToby Isaac   ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr);
1509f1a1a6beSBarry Smith   ierr = PetscOptionsBoundedInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg,0);CHKERRQ(ierr);
1510dd8e54a2SToby Isaac   if (flg) {
1511dd8e54a2SToby Isaac     ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr);
1512f885a11aSToby Isaac   } else {
1513dd8e54a2SToby Isaac     ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr);
15145a856986SBarry Smith     ierr = PetscOptionsBoundedInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg,1);CHKERRQ(ierr);
1515dd8e54a2SToby Isaac     if (flg) {
1516dd8e54a2SToby Isaac       ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr);
1517dd8e54a2SToby Isaac     }
1518dd8e54a2SToby Isaac   }
1519dd8e54a2SToby Isaac   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
15205a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg,0);CHKERRQ(ierr);
1521dd8e54a2SToby Isaac   if (flg) {
1522dd8e54a2SToby Isaac     ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr);
1523dd8e54a2SToby Isaac   }
1524a6121fbdSMatthew G. Knepley #if 0
15255a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-dm_refine","equivalent to -dm_forest_set_minimum_refinement and -dm_forest_set_initial_refinement with the same value",NULL,minRefinement,&minRefinement,&flg,0);CHKERRQ(ierr);
1526a6121fbdSMatthew G. Knepley   if (flg) {
1527a6121fbdSMatthew G. Knepley     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1528a6121fbdSMatthew G. Knepley     ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr);
1529a6121fbdSMatthew G. Knepley   }
15305a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-dm_refine_hierarchy","equivalent to -dm_forest_set_minimum_refinement 0 and -dm_forest_set_initial_refinement",NULL,initRefinement,&initRefinement,&flg,0);CHKERRQ(ierr);
1531a6121fbdSMatthew G. Knepley   if (flg) {
1532a6121fbdSMatthew G. Knepley     ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr);
1533a6121fbdSMatthew G. Knepley     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
1534a6121fbdSMatthew G. Knepley   }
1535a6121fbdSMatthew G. Knepley #endif
1536dd8e54a2SToby Isaac   ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr);
15375a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg,0);CHKERRQ(ierr);
1538dd8e54a2SToby Isaac   if (flg) {
1539dd8e54a2SToby Isaac     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1540db4d5e8cSToby Isaac   }
154156ba9f64SToby Isaac   ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr);
15425a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg,0);CHKERRQ(ierr);
154356ba9f64SToby Isaac   if (flg) {
154456ba9f64SToby Isaac     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
154556ba9f64SToby Isaac   }
1546c7eeac06SToby Isaac   ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr);
15475a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg,0);CHKERRQ(ierr);
1548c7eeac06SToby Isaac   if (flg) {
1549c7eeac06SToby Isaac     ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr);
1550c7eeac06SToby Isaac   }
1551c7eeac06SToby Isaac   ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr);
1552589a23caSBarry Smith   ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,sizeof(stringBuffer),&flg);CHKERRQ(ierr);
1553c7eeac06SToby Isaac   if (flg) {
1554c7eeac06SToby Isaac     ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr);
1555c7eeac06SToby Isaac   }
1556c7eeac06SToby Isaac   ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr);
15575a856986SBarry Smith   ierr = PetscOptionsBoundedInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg,0);CHKERRQ(ierr);
1558c7eeac06SToby Isaac   if (flg) {
1559c7eeac06SToby Isaac     ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr);
1560c7eeac06SToby Isaac   }
1561c7eeac06SToby Isaac   ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr);
1562c7eeac06SToby Isaac   ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr);
1563c7eeac06SToby Isaac   if (flg) {
1564c7eeac06SToby Isaac     ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr);
1565c7eeac06SToby Isaac   }
1566db4d5e8cSToby Isaac   ierr = PetscOptionsTail();CHKERRQ(ierr);
1567db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1568db4d5e8cSToby Isaac }
1569db4d5e8cSToby Isaac 
1570276c5506SMatthew G. Knepley PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1571d8984e3bSMatthew G. Knepley {
1572d8984e3bSMatthew G. Knepley   PetscErrorCode ierr;
1573d8984e3bSMatthew G. Knepley 
1574d8984e3bSMatthew G. Knepley   PetscFunctionBegin;
1575d8984e3bSMatthew G. Knepley   if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);}
1576792b654fSMatthew G. Knepley   ierr = DMCreateSectionSubDM(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
1577d8984e3bSMatthew G. Knepley   PetscFunctionReturn(0);
1578d8984e3bSMatthew G. Knepley }
1579d8984e3bSMatthew G. Knepley 
15805421bac9SToby Isaac PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
15815421bac9SToby Isaac {
15825421bac9SToby Isaac   DMLabel        refine;
15835421bac9SToby Isaac   DM             fineDM;
15845421bac9SToby Isaac   PetscErrorCode ierr;
15855421bac9SToby Isaac 
15865421bac9SToby Isaac   PetscFunctionBegin;
15875421bac9SToby Isaac   ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr);
15885421bac9SToby Isaac   if (fineDM) {
15895421bac9SToby Isaac     ierr       = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr);
15905421bac9SToby Isaac     *dmRefined = fineDM;
15915421bac9SToby Isaac     PetscFunctionReturn(0);
15925421bac9SToby Isaac   }
15935421bac9SToby Isaac   ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr);
15945421bac9SToby Isaac   ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr);
15955421bac9SToby Isaac   if (!refine) {
1596d67d17b1SMatthew G. Knepley     ierr = DMLabelCreate(PETSC_COMM_SELF, "refine",&refine);CHKERRQ(ierr);
1597a1b0c543SToby Isaac     ierr = DMLabelSetDefaultValue(refine,DM_ADAPT_REFINE);CHKERRQ(ierr);
1598d67d17b1SMatthew G. Knepley   } else {
1599d67d17b1SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) refine);CHKERRQ(ierr);
1600a1b0c543SToby Isaac   }
1601a1b0c543SToby Isaac   ierr = DMForestSetAdaptivityLabel(*dmRefined,refine);CHKERRQ(ierr);
1602a1b0c543SToby Isaac   ierr = DMLabelDestroy(&refine);CHKERRQ(ierr);
16035421bac9SToby Isaac   PetscFunctionReturn(0);
16045421bac9SToby Isaac }
16055421bac9SToby Isaac 
16065421bac9SToby Isaac PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
16075421bac9SToby Isaac {
16085421bac9SToby Isaac   DMLabel        coarsen;
16095421bac9SToby Isaac   DM             coarseDM;
16105421bac9SToby Isaac   PetscErrorCode ierr;
16115421bac9SToby Isaac 
16125421bac9SToby Isaac   PetscFunctionBegin;
16134098eed7SToby Isaac   {
16144098eed7SToby Isaac     PetscMPIInt mpiComparison;
16154098eed7SToby Isaac     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);
16164098eed7SToby Isaac 
1617ffc4695bSBarry Smith     ierr = MPI_Comm_compare(comm, dmcomm, &mpiComparison);CHKERRMPI(ierr);
1618f885a11aSToby Isaac     if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet");
16194098eed7SToby Isaac   }
16205421bac9SToby Isaac   ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr);
16215421bac9SToby Isaac   if (coarseDM) {
16225421bac9SToby Isaac     ierr         = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr);
16235421bac9SToby Isaac     *dmCoarsened = coarseDM;
16245421bac9SToby Isaac     PetscFunctionReturn(0);
16255421bac9SToby Isaac   }
16265421bac9SToby Isaac   ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr);
16277ee90a76SStefano Zampini   ierr = DMForestSetAdaptivityPurpose(*dmCoarsened,DM_ADAPT_COARSEN);CHKERRQ(ierr);
16285421bac9SToby Isaac   ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr);
16295421bac9SToby Isaac   if (!coarsen) {
1630d67d17b1SMatthew G. Knepley     ierr = DMLabelCreate(PETSC_COMM_SELF, "coarsen",&coarsen);CHKERRQ(ierr);
1631a1b0c543SToby Isaac     ierr = DMLabelSetDefaultValue(coarsen,DM_ADAPT_COARSEN);CHKERRQ(ierr);
1632a1b0c543SToby Isaac   } else {
1633d67d17b1SMatthew G. Knepley     ierr = PetscObjectReference((PetscObject) coarsen);CHKERRQ(ierr);
16345421bac9SToby Isaac   }
1635a1b0c543SToby Isaac   ierr = DMForestSetAdaptivityLabel(*dmCoarsened,coarsen);CHKERRQ(ierr);
1636a1b0c543SToby Isaac   ierr = DMLabelDestroy(&coarsen);CHKERRQ(ierr);
16375421bac9SToby Isaac   PetscFunctionReturn(0);
16385421bac9SToby Isaac }
16395421bac9SToby Isaac 
1640a1b0c543SToby Isaac static PetscErrorCode DMAdaptLabel_Forest(DM dm, DMLabel label, DM *adaptedDM)
164109350103SToby Isaac {
164209350103SToby Isaac   PetscBool      success;
164309350103SToby Isaac   PetscErrorCode ierr;
164409350103SToby Isaac 
164509350103SToby Isaac   PetscFunctionBegin;
164609350103SToby Isaac   ierr = DMForestTemplate(dm,PetscObjectComm((PetscObject)dm),adaptedDM);CHKERRQ(ierr);
1647a1b0c543SToby Isaac   ierr = DMForestSetAdaptivityLabel(*adaptedDM,label);CHKERRQ(ierr);
164809350103SToby Isaac   ierr = DMSetUp(*adaptedDM);CHKERRQ(ierr);
164909350103SToby Isaac   ierr = DMForestGetAdaptivitySuccess(*adaptedDM,&success);CHKERRQ(ierr);
165009350103SToby Isaac   if (!success) {
165109350103SToby Isaac     ierr = DMDestroy(adaptedDM);CHKERRQ(ierr);
165209350103SToby Isaac     *adaptedDM = NULL;
165309350103SToby Isaac   }
165409350103SToby Isaac   PetscFunctionReturn(0);
165509350103SToby Isaac }
165609350103SToby Isaac 
1657d222f98bSToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm)
1658d222f98bSToby Isaac {
1659d222f98bSToby Isaac   PetscErrorCode ierr;
1660d222f98bSToby Isaac 
1661d222f98bSToby Isaac   PetscFunctionBegin;
1662d222f98bSToby Isaac   ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr);
1663d222f98bSToby Isaac 
1664d222f98bSToby Isaac   dm->ops->clone          = DMClone_Forest;
1665d222f98bSToby Isaac   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1666d222f98bSToby Isaac   dm->ops->destroy        = DMDestroy_Forest;
1667d8984e3bSMatthew G. Knepley   dm->ops->createsubdm    = DMCreateSubDM_Forest;
16685421bac9SToby Isaac   dm->ops->refine         = DMRefine_Forest;
16695421bac9SToby Isaac   dm->ops->coarsen        = DMCoarsen_Forest;
16700d1cd5e0SMatthew G. Knepley   dm->ops->adaptlabel     = DMAdaptLabel_Forest;
1671d222f98bSToby Isaac   PetscFunctionReturn(0);
1672d222f98bSToby Isaac }
1673d222f98bSToby Isaac 
16749be51f97SToby Isaac /*MC
16759be51f97SToby Isaac 
1676bae1f979SBarry Smith      DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh.  Forests usually have a base DM
1677bae1f979SBarry Smith   (see DMForestGetBaseDM()), from which it is refined.  The refinement and partitioning of forests is considered
1678bae1f979SBarry Smith   immutable after DMSetUp() is called.  To adapt a mesh, one should call DMForestTemplate() to create a new mesh that
1679bae1f979SBarry Smith   will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new
1680bae1f979SBarry Smith   mesh.
1681bae1f979SBarry Smith 
1682bae1f979SBarry Smith   To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the
1683bae1f979SBarry Smith   previous mesh whose values indicate which cells should be refined (DM_ADAPT_REFINE) or coarsened (DM_ADAPT_COARSEN)
1684bae1f979SBarry Smith   and how (subtypes are free to allow additional values for things like anisotropic refinement).  The label should be
1685bae1f979SBarry Smith   given to the *new* mesh with DMForestSetAdaptivityLabel().
16869be51f97SToby Isaac 
16879be51f97SToby Isaac   Level: advanced
16889be51f97SToby Isaac 
16899be51f97SToby Isaac .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel()
16909be51f97SToby Isaac M*/
16919be51f97SToby Isaac 
1692db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1693db4d5e8cSToby Isaac {
1694db4d5e8cSToby Isaac   DM_Forest      *forest;
1695db4d5e8cSToby Isaac   PetscErrorCode ierr;
1696db4d5e8cSToby Isaac 
1697db4d5e8cSToby Isaac   PetscFunctionBegin;
1698db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1699db4d5e8cSToby Isaac   ierr                         = PetscNewLog(dm,&forest);CHKERRQ(ierr);
1700db4d5e8cSToby Isaac   dm->dim                      = 0;
1701db4d5e8cSToby Isaac   dm->data                     = forest;
1702db4d5e8cSToby Isaac   forest->refct                = 1;
1703db4d5e8cSToby Isaac   forest->data                 = NULL;
1704db4d5e8cSToby Isaac   forest->topology             = NULL;
1705cd3c525cSToby Isaac   forest->adapt                = NULL;
1706db4d5e8cSToby Isaac   forest->base                 = NULL;
17076a87ffbfSToby Isaac   forest->adaptPurpose         = DM_ADAPT_DETERMINE;
1708db4d5e8cSToby Isaac   forest->adjDim               = PETSC_DEFAULT;
1709db4d5e8cSToby Isaac   forest->overlap              = PETSC_DEFAULT;
1710db4d5e8cSToby Isaac   forest->minRefinement        = PETSC_DEFAULT;
1711db4d5e8cSToby Isaac   forest->maxRefinement        = PETSC_DEFAULT;
171256ba9f64SToby Isaac   forest->initRefinement       = PETSC_DEFAULT;
1713c7eeac06SToby Isaac   forest->cStart               = PETSC_DETERMINE;
1714c7eeac06SToby Isaac   forest->cEnd                 = PETSC_DETERMINE;
1715cd3c525cSToby Isaac   forest->cellSF               = NULL;
1716ebdf65a2SToby Isaac   forest->adaptLabel           = NULL;
1717db4d5e8cSToby Isaac   forest->gradeFactor          = 2;
1718db4d5e8cSToby Isaac   forest->cellWeights          = NULL;
1719db4d5e8cSToby Isaac   forest->cellWeightsCopyMode  = PETSC_USE_POINTER;
1720db4d5e8cSToby Isaac   forest->weightsFactor        = 1.;
1721db4d5e8cSToby Isaac   forest->weightCapacity       = 1.;
1722a73e2921SToby Isaac   ierr                         = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr);
1723d222f98bSToby Isaac   ierr                         = DMInitialize_Forest(dm);CHKERRQ(ierr);
1724db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1725db4d5e8cSToby Isaac }
1726