xref: /petsc/src/dm/impls/plex/plexrefine.c (revision dadcf80911fb48939c55327431ae8d7e47dbe367)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2412e9a14SMatthew G. Knepley #include <petsc/private/petscfeimpl.h>  /* For PetscFEInterpolate_Static() */
3012bc364SMatthew G. Knepley 
4012bc364SMatthew G. Knepley #include <petscdmplextransform.h>
575d3a19aSMatthew G. Knepley #include <petscsf.h>
675d3a19aSMatthew G. Knepley 
7963fc26aSMatthew G. Knepley /*@
8963fc26aSMatthew G. Knepley   DMPlexCreateProcessSF - Create an SF which just has process connectivity
9963fc26aSMatthew G. Knepley 
10d083f849SBarry Smith   Collective on dm
11963fc26aSMatthew G. Knepley 
12963fc26aSMatthew G. Knepley   Input Parameters:
13963fc26aSMatthew G. Knepley + dm      - The DM
14963fc26aSMatthew G. Knepley - sfPoint - The PetscSF which encodes point connectivity
15963fc26aSMatthew G. Knepley 
16963fc26aSMatthew G. Knepley   Output Parameters:
17963fc26aSMatthew G. Knepley + processRanks - A list of process neighbors, or NULL
18963fc26aSMatthew G. Knepley - sfProcess    - An SF encoding the process connectivity, or NULL
19963fc26aSMatthew G. Knepley 
20963fc26aSMatthew G. Knepley   Level: developer
21963fc26aSMatthew G. Knepley 
22963fc26aSMatthew G. Knepley .seealso: PetscSFCreate(), DMPlexCreateTwoSidedProcessSF()
23963fc26aSMatthew G. Knepley @*/
24963fc26aSMatthew G. Knepley PetscErrorCode DMPlexCreateProcessSF(DM dm, PetscSF sfPoint, IS *processRanks, PetscSF *sfProcess)
2575d3a19aSMatthew G. Knepley {
2675d3a19aSMatthew G. Knepley   PetscInt           numRoots, numLeaves, l;
2775d3a19aSMatthew G. Knepley   const PetscInt    *localPoints;
2875d3a19aSMatthew G. Knepley   const PetscSFNode *remotePoints;
2975d3a19aSMatthew G. Knepley   PetscInt          *localPointsNew;
3075d3a19aSMatthew G. Knepley   PetscSFNode       *remotePointsNew;
3175d3a19aSMatthew G. Knepley   PetscInt          *ranks, *ranksNew;
329852e123SBarry Smith   PetscMPIInt        size;
3375d3a19aSMatthew G. Knepley 
3475d3a19aSMatthew G. Knepley   PetscFunctionBegin;
35963fc26aSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
36963fc26aSMatthew G. Knepley   PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2);
37963fc26aSMatthew G. Knepley   if (processRanks) {PetscValidPointer(processRanks, 3);}
38963fc26aSMatthew G. Knepley   if (sfProcess)    {PetscValidPointer(sfProcess, 4);}
395f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(numLeaves, &ranks));
4275d3a19aSMatthew G. Knepley   for (l = 0; l < numLeaves; ++l) {
4375d3a19aSMatthew G. Knepley     ranks[l] = remotePoints[l].rank;
4475d3a19aSMatthew G. Knepley   }
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSortRemoveDupsInt(&numLeaves, ranks));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(numLeaves, &ranksNew));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(numLeaves, &localPointsNew));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(numLeaves, &remotePointsNew));
4975d3a19aSMatthew G. Knepley   for (l = 0; l < numLeaves; ++l) {
5075d3a19aSMatthew G. Knepley     ranksNew[l]              = ranks[l];
5175d3a19aSMatthew G. Knepley     localPointsNew[l]        = l;
5275d3a19aSMatthew G. Knepley     remotePointsNew[l].index = 0;
5375d3a19aSMatthew G. Knepley     remotePointsNew[l].rank  = ranksNew[l];
5475d3a19aSMatthew G. Knepley   }
555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(ranks));
565f80ce2aSJacob Faibussowitsch   if (processRanks) CHKERRQ(ISCreateGeneral(PetscObjectComm((PetscObject)dm), numLeaves, ranksNew, PETSC_OWN_POINTER, processRanks));
575f80ce2aSJacob Faibussowitsch   else              CHKERRQ(PetscFree(ranksNew));
58963fc26aSMatthew G. Knepley   if (sfProcess) {
595f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess));
605f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetName((PetscObject) *sfProcess, "Process SF"));
615f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFSetFromOptions(*sfProcess));
625f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFSetGraph(*sfProcess, size, numLeaves, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER));
63963fc26aSMatthew G. Knepley   }
6475d3a19aSMatthew G. Knepley   PetscFunctionReturn(0);
6575d3a19aSMatthew G. Knepley }
6675d3a19aSMatthew G. Knepley 
672389894bSMatthew G. Knepley /*@
682389894bSMatthew G. Knepley   DMPlexCreateCoarsePointIS - Creates an IS covering the coarse DM chart with the fine points as data
692389894bSMatthew G. Knepley 
70012bc364SMatthew G. Knepley   Collective on dm
71012bc364SMatthew G. Knepley 
722389894bSMatthew G. Knepley   Input Parameter:
732389894bSMatthew G. Knepley . dm - The coarse DM
742389894bSMatthew G. Knepley 
752389894bSMatthew G. Knepley   Output Parameter:
762389894bSMatthew G. Knepley . fpointIS - The IS of all the fine points which exist in the original coarse mesh
772389894bSMatthew G. Knepley 
782389894bSMatthew G. Knepley   Level: developer
792389894bSMatthew G. Knepley 
8097d8846cSMatthew Knepley .seealso: DMRefine(), DMPlexSetRefinementUniform(), DMPlexGetSubpointIS()
812389894bSMatthew G. Knepley @*/
822389894bSMatthew G. Knepley PetscErrorCode DMPlexCreateCoarsePointIS(DM dm, IS *fpointIS)
832389894bSMatthew G. Knepley {
84012bc364SMatthew G. Knepley   DMPlexTransform tr;
85412e9a14SMatthew G. Knepley   PetscInt       *fpoints;
86327c2912SStefano Zampini   PetscInt        pStart, pEnd, p, vStart, vEnd, v;
872389894bSMatthew G. Knepley 
882389894bSMatthew G. Knepley   PetscFunctionBegin;
895f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetChart(dm, &pStart, &pEnd));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexTransformCreate(PetscObjectComm((PetscObject) dm), &tr));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexTransformSetUp(tr));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(pEnd-pStart, &fpoints));
942389894bSMatthew G. Knepley   for (p = 0; p < pEnd-pStart; ++p) fpoints[p] = -1;
95412e9a14SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
96012bc364SMatthew G. Knepley     PetscInt vNew = -1; /* quiet overzealous may be used uninitialized check */
97327c2912SStefano Zampini 
985f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexTransformGetTargetPoint(tr, DM_POLYTOPE_POINT, DM_POLYTOPE_POINT, p, 0, &vNew));
99412e9a14SMatthew G. Knepley     fpoints[v-pStart] = vNew;
1002389894bSMatthew G. Knepley   }
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexTransformDestroy(&tr));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF, pEnd-pStart, fpoints, PETSC_OWN_POINTER, fpointIS));
1032389894bSMatthew G. Knepley   PetscFunctionReturn(0);
1042389894bSMatthew G. Knepley }
1052389894bSMatthew G. Knepley 
106012bc364SMatthew G. Knepley /*@C
107012bc364SMatthew G. Knepley   DMPlexSetTransformType - Set the transform type for uniform refinement
108012bc364SMatthew G. Knepley 
109012bc364SMatthew G. Knepley   Input Parameters:
110012bc364SMatthew G. Knepley + dm - The DM
111012bc364SMatthew G. Knepley - type - The transform type for uniform refinement
112012bc364SMatthew G. Knepley 
113012bc364SMatthew G. Knepley   Level: developer
114012bc364SMatthew G. Knepley 
115012bc364SMatthew G. Knepley .seealso: DMPlexTransformType, DMRefine(), DMPlexGetTransformType(), DMPlexSetRefinementUniform()
116012bc364SMatthew G. Knepley @*/
117012bc364SMatthew G. Knepley PetscErrorCode DMPlexSetTransformType(DM dm, DMPlexTransformType type)
118012bc364SMatthew G. Knepley {
119012bc364SMatthew G. Knepley   DM_Plex        *mesh = (DM_Plex*) dm->data;
120012bc364SMatthew G. Knepley 
121012bc364SMatthew G. Knepley   PetscFunctionBegin;
122012bc364SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
123012bc364SMatthew G. Knepley   if (type) PetscValidCharPointer(type, 2);
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(mesh->transformType));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrallocpy(type, &mesh->transformType));
126012bc364SMatthew G. Knepley   PetscFunctionReturn(0);
127012bc364SMatthew G. Knepley }
128012bc364SMatthew G. Knepley 
129012bc364SMatthew G. Knepley /*@C
130012bc364SMatthew G. Knepley   DMPlexGetTransformType - Retrieve the transform type for uniform refinement
131012bc364SMatthew G. Knepley 
132012bc364SMatthew G. Knepley   Input Parameter:
133012bc364SMatthew G. Knepley . dm - The DM
134012bc364SMatthew G. Knepley 
135012bc364SMatthew G. Knepley   Output Parameter:
136012bc364SMatthew G. Knepley . type - The transform type for uniform refinement
137012bc364SMatthew G. Knepley 
138012bc364SMatthew G. Knepley   Level: developer
139012bc364SMatthew G. Knepley 
140012bc364SMatthew G. Knepley .seealso: DMPlexTransformType, DMRefine(), DMPlexSetTransformType(), DMPlexGetRefinementUniform()
141012bc364SMatthew G. Knepley @*/
142012bc364SMatthew G. Knepley PetscErrorCode DMPlexGetTransformType(DM dm, DMPlexTransformType *type)
143012bc364SMatthew G. Knepley {
144012bc364SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
145012bc364SMatthew G. Knepley 
146012bc364SMatthew G. Knepley   PetscFunctionBegin;
147012bc364SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
148012bc364SMatthew G. Knepley   PetscValidPointer(type, 2);
149012bc364SMatthew G. Knepley   *type = mesh->transformType;
150012bc364SMatthew G. Knepley   PetscFunctionReturn(0);
151012bc364SMatthew G. Knepley }
152012bc364SMatthew G. Knepley 
1530e2b6761SMatthew G. Knepley /*@
1540e2b6761SMatthew G. Knepley   DMPlexSetRefinementUniform - Set the flag for uniform refinement
1550e2b6761SMatthew G. Knepley 
1560e2b6761SMatthew G. Knepley   Input Parameters:
1570e2b6761SMatthew G. Knepley + dm - The DM
1580e2b6761SMatthew G. Knepley - refinementUniform - The flag for uniform refinement
1590e2b6761SMatthew G. Knepley 
1600e2b6761SMatthew G. Knepley   Level: developer
1610e2b6761SMatthew G. Knepley 
1620e2b6761SMatthew G. Knepley .seealso: DMRefine(), DMPlexGetRefinementUniform(), DMPlexGetRefinementLimit(), DMPlexSetRefinementLimit()
1630e2b6761SMatthew G. Knepley @*/
16475d3a19aSMatthew G. Knepley PetscErrorCode DMPlexSetRefinementUniform(DM dm, PetscBool refinementUniform)
16575d3a19aSMatthew G. Knepley {
16675d3a19aSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
16775d3a19aSMatthew G. Knepley 
16875d3a19aSMatthew G. Knepley   PetscFunctionBegin;
169012bc364SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
17075d3a19aSMatthew G. Knepley   mesh->refinementUniform = refinementUniform;
17175d3a19aSMatthew G. Knepley   PetscFunctionReturn(0);
17275d3a19aSMatthew G. Knepley }
17375d3a19aSMatthew G. Knepley 
1740e2b6761SMatthew G. Knepley /*@
1750e2b6761SMatthew G. Knepley   DMPlexGetRefinementUniform - Retrieve the flag for uniform refinement
1760e2b6761SMatthew G. Knepley 
1770e2b6761SMatthew G. Knepley   Input Parameter:
1780e2b6761SMatthew G. Knepley . dm - The DM
1790e2b6761SMatthew G. Knepley 
1800e2b6761SMatthew G. Knepley   Output Parameter:
1810e2b6761SMatthew G. Knepley . refinementUniform - The flag for uniform refinement
1820e2b6761SMatthew G. Knepley 
1830e2b6761SMatthew G. Knepley   Level: developer
1840e2b6761SMatthew G. Knepley 
1850e2b6761SMatthew G. Knepley .seealso: DMRefine(), DMPlexSetRefinementUniform(), DMPlexGetRefinementLimit(), DMPlexSetRefinementLimit()
1860e2b6761SMatthew G. Knepley @*/
18775d3a19aSMatthew G. Knepley PetscErrorCode DMPlexGetRefinementUniform(DM dm, PetscBool *refinementUniform)
18875d3a19aSMatthew G. Knepley {
18975d3a19aSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
19075d3a19aSMatthew G. Knepley 
19175d3a19aSMatthew G. Knepley   PetscFunctionBegin;
192012bc364SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
193*dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(refinementUniform,  2);
19475d3a19aSMatthew G. Knepley   *refinementUniform = mesh->refinementUniform;
19575d3a19aSMatthew G. Knepley   PetscFunctionReturn(0);
19675d3a19aSMatthew G. Knepley }
19775d3a19aSMatthew G. Knepley 
1980e2b6761SMatthew G. Knepley /*@
1990e2b6761SMatthew G. Knepley   DMPlexSetRefinementLimit - Set the maximum cell volume for refinement
2000e2b6761SMatthew G. Knepley 
2010e2b6761SMatthew G. Knepley   Input Parameters:
2020e2b6761SMatthew G. Knepley + dm - The DM
2030e2b6761SMatthew G. Knepley - refinementLimit - The maximum cell volume in the refined mesh
2040e2b6761SMatthew G. Knepley 
2050e2b6761SMatthew G. Knepley   Level: developer
2060e2b6761SMatthew G. Knepley 
2070e2b6761SMatthew G. Knepley .seealso: DMRefine(), DMPlexGetRefinementLimit(), DMPlexGetRefinementUniform(), DMPlexSetRefinementUniform()
2080e2b6761SMatthew G. Knepley @*/
20975d3a19aSMatthew G. Knepley PetscErrorCode DMPlexSetRefinementLimit(DM dm, PetscReal refinementLimit)
21075d3a19aSMatthew G. Knepley {
21175d3a19aSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
21275d3a19aSMatthew G. Knepley 
21375d3a19aSMatthew G. Knepley   PetscFunctionBegin;
21475d3a19aSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
21575d3a19aSMatthew G. Knepley   mesh->refinementLimit = refinementLimit;
21675d3a19aSMatthew G. Knepley   PetscFunctionReturn(0);
21775d3a19aSMatthew G. Knepley }
21875d3a19aSMatthew G. Knepley 
2190e2b6761SMatthew G. Knepley /*@
2200e2b6761SMatthew G. Knepley   DMPlexGetRefinementLimit - Retrieve the maximum cell volume for refinement
2210e2b6761SMatthew G. Knepley 
2220e2b6761SMatthew G. Knepley   Input Parameter:
2230e2b6761SMatthew G. Knepley . dm - The DM
2240e2b6761SMatthew G. Knepley 
2250e2b6761SMatthew G. Knepley   Output Parameter:
2260e2b6761SMatthew G. Knepley . refinementLimit - The maximum cell volume in the refined mesh
2270e2b6761SMatthew G. Knepley 
2280e2b6761SMatthew G. Knepley   Level: developer
2290e2b6761SMatthew G. Knepley 
2300e2b6761SMatthew G. Knepley .seealso: DMRefine(), DMPlexSetRefinementLimit(), DMPlexGetRefinementUniform(), DMPlexSetRefinementUniform()
2310e2b6761SMatthew G. Knepley @*/
23275d3a19aSMatthew G. Knepley PetscErrorCode DMPlexGetRefinementLimit(DM dm, PetscReal *refinementLimit)
23375d3a19aSMatthew G. Knepley {
23475d3a19aSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
23575d3a19aSMatthew G. Knepley 
23675d3a19aSMatthew G. Knepley   PetscFunctionBegin;
23775d3a19aSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
238*dadcf809SJacob Faibussowitsch   PetscValidRealPointer(refinementLimit,  2);
23975d3a19aSMatthew G. Knepley   /* if (mesh->refinementLimit < 0) = getMaxVolume()/2.0; */
24075d3a19aSMatthew G. Knepley   *refinementLimit = mesh->refinementLimit;
24175d3a19aSMatthew G. Knepley   PetscFunctionReturn(0);
24275d3a19aSMatthew G. Knepley }
24375d3a19aSMatthew G. Knepley 
244b28003e6SMatthew G. Knepley /*@
245b28003e6SMatthew G. Knepley   DMPlexSetRefinementFunction - Set the function giving the maximum cell volume for refinement
246b28003e6SMatthew G. Knepley 
247b28003e6SMatthew G. Knepley   Input Parameters:
248b28003e6SMatthew G. Knepley + dm - The DM
249b28003e6SMatthew G. Knepley - refinementFunc - Function giving the maximum cell volume in the refined mesh
250b28003e6SMatthew G. Knepley 
251b28003e6SMatthew G. Knepley   Note: The calling sequence is refinementFunc(coords, limit)
252b28003e6SMatthew G. Knepley $ coords - Coordinates of the current point, usually a cell centroid
253b28003e6SMatthew G. Knepley $ limit  - The maximum cell volume for a cell containing this point
254b28003e6SMatthew G. Knepley 
255b28003e6SMatthew G. Knepley   Level: developer
256b28003e6SMatthew G. Knepley 
257b28003e6SMatthew G. Knepley .seealso: DMRefine(), DMPlexGetRefinementFunction(), DMPlexGetRefinementUniform(), DMPlexSetRefinementUniform(), DMPlexGetRefinementLimit(), DMPlexSetRefinementLimit()
258b28003e6SMatthew G. Knepley @*/
259b28003e6SMatthew G. Knepley PetscErrorCode DMPlexSetRefinementFunction(DM dm, PetscErrorCode (*refinementFunc)(const PetscReal [], PetscReal *))
260b28003e6SMatthew G. Knepley {
261b28003e6SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
262b28003e6SMatthew G. Knepley 
263b28003e6SMatthew G. Knepley   PetscFunctionBegin;
264b28003e6SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
265b28003e6SMatthew G. Knepley   mesh->refinementFunc = refinementFunc;
266b28003e6SMatthew G. Knepley   PetscFunctionReturn(0);
267b28003e6SMatthew G. Knepley }
268b28003e6SMatthew G. Knepley 
269b28003e6SMatthew G. Knepley /*@
270b28003e6SMatthew G. Knepley   DMPlexGetRefinementFunction - Get the function giving the maximum cell volume for refinement
271b28003e6SMatthew G. Knepley 
272b28003e6SMatthew G. Knepley   Input Parameter:
273b28003e6SMatthew G. Knepley . dm - The DM
274b28003e6SMatthew G. Knepley 
275b28003e6SMatthew G. Knepley   Output Parameter:
276b28003e6SMatthew G. Knepley . refinementFunc - Function giving the maximum cell volume in the refined mesh
277b28003e6SMatthew G. Knepley 
278b28003e6SMatthew G. Knepley   Note: The calling sequence is refinementFunc(coords, limit)
279b28003e6SMatthew G. Knepley $ coords - Coordinates of the current point, usually a cell centroid
280b28003e6SMatthew G. Knepley $ limit  - The maximum cell volume for a cell containing this point
281b28003e6SMatthew G. Knepley 
282b28003e6SMatthew G. Knepley   Level: developer
283b28003e6SMatthew G. Knepley 
284b28003e6SMatthew G. Knepley .seealso: DMRefine(), DMPlexSetRefinementFunction(), DMPlexGetRefinementUniform(), DMPlexSetRefinementUniform(), DMPlexGetRefinementLimit(), DMPlexSetRefinementLimit()
285b28003e6SMatthew G. Knepley @*/
286b28003e6SMatthew G. Knepley PetscErrorCode DMPlexGetRefinementFunction(DM dm, PetscErrorCode (**refinementFunc)(const PetscReal [], PetscReal *))
287b28003e6SMatthew G. Knepley {
288b28003e6SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
289b28003e6SMatthew G. Knepley 
290b28003e6SMatthew G. Knepley   PetscFunctionBegin;
291b28003e6SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
292b28003e6SMatthew G. Knepley   PetscValidPointer(refinementFunc,  2);
293b28003e6SMatthew G. Knepley   *refinementFunc = mesh->refinementFunc;
294b28003e6SMatthew G. Knepley   PetscFunctionReturn(0);
295b28003e6SMatthew G. Knepley }
296b28003e6SMatthew G. Knepley 
297012bc364SMatthew G. Knepley PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *rdm)
2980d1cd5e0SMatthew G. Knepley {
299492b8470SStefano Zampini   PetscBool      isUniform;
3000d1cd5e0SMatthew G. Knepley 
3010d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
3025f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetRefinementUniform(dm, &isUniform));
3035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm, NULL, "-initref_dm_view"));
3040d1cd5e0SMatthew G. Knepley   if (isUniform) {
305012bc364SMatthew G. Knepley     DMPlexTransform     tr;
30651a74b61SMatthew G. Knepley     DM                  cdm, rcdm;
307012bc364SMatthew G. Knepley     DMPlexTransformType trType;
308012bc364SMatthew G. Knepley     const char         *prefix;
309012bc364SMatthew G. Knepley     PetscOptions       options;
3100d1cd5e0SMatthew G. Knepley 
3115f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexTransformCreate(PetscObjectComm((PetscObject) dm), &tr));
3125f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexTransformSetDM(tr, dm));
3135f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetTransformType(dm, &trType));
3145f80ce2aSJacob Faibussowitsch     if (trType) CHKERRQ(DMPlexTransformSetType(tr, trType));
3155f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject) dm, &prefix));
3165f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) tr,  prefix));
3175f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectGetOptions((PetscObject) dm, &options));
3185f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetOptions((PetscObject) tr, options));
3195f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexTransformSetFromOptions(tr));
3205f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectSetOptions((PetscObject) tr, NULL));
3215f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexTransformSetUp(tr));
3225f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectViewFromOptions((PetscObject) tr, NULL, "-dm_plex_transform_view"));
3235f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexTransformApply(tr, dm, rdm));
3245f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexSetRegularRefinement(*rdm, PETSC_TRUE));
3255f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCopyDisc(dm, *rdm));
3265f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinateDM(dm, &cdm));
3275f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetCoordinateDM(*rdm, &rcdm));
3285f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCopyDisc(cdm, rcdm));
3295f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexTransformCreateDiscLabels(tr, *rdm));
3305f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexTransformDestroy(&tr));
3310d1cd5e0SMatthew G. Knepley   } else {
3325f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexRefine_Internal(dm, NULL, NULL, NULL, rdm));
3330d1cd5e0SMatthew G. Knepley   }
334012bc364SMatthew G. Knepley   if (*rdm) {
335012bc364SMatthew G. Knepley     ((DM_Plex *) (*rdm)->data)->printFEM = ((DM_Plex *) dm->data)->printFEM;
336012bc364SMatthew G. Knepley     ((DM_Plex *) (*rdm)->data)->printL2  = ((DM_Plex *) dm->data)->printL2;
3376d7c9049SMatthew G. Knepley   }
3385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm, NULL, "-postref_dm_view"));
3390d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
3400d1cd5e0SMatthew G. Knepley }
3410d1cd5e0SMatthew G. Knepley 
342012bc364SMatthew G. Knepley PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM rdm[])
3430d1cd5e0SMatthew G. Knepley {
3440d1cd5e0SMatthew G. Knepley   DM             cdm = dm;
3450d1cd5e0SMatthew G. Knepley   PetscInt       r;
3460d1cd5e0SMatthew G. Knepley   PetscBool      isUniform, localized;
3470d1cd5e0SMatthew G. Knepley 
3480d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
3495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetRefinementUniform(dm, &isUniform));
3505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocalized(dm, &localized));
3510d1cd5e0SMatthew G. Knepley   if (isUniform) {
3520d1cd5e0SMatthew G. Knepley     for (r = 0; r < nlevels; ++r) {
353012bc364SMatthew G. Knepley       DMPlexTransform tr;
35451a74b61SMatthew G. Knepley       DM              codm, rcodm;
355012bc364SMatthew G. Knepley       const char     *prefix;
3560d1cd5e0SMatthew G. Knepley 
3575f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexTransformCreate(PetscObjectComm((PetscObject) cdm), &tr));
3585f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject) cdm, &prefix));
3595f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject) tr,   prefix));
3605f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexTransformSetDM(tr, cdm));
3615f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexTransformSetFromOptions(tr));
3625f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexTransformSetUp(tr));
3635f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexTransformApply(tr, cdm, &rdm[r]));
3645f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetCoarsenLevel(rdm[r], cdm->leveldown));
3655f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetRefineLevel(rdm[r], cdm->levelup+1));
3665f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCopyDisc(cdm, rdm[r]));
3675f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetCoordinateDM(dm, &codm));
3685f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetCoordinateDM(rdm[r], &rcodm));
3695f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCopyDisc(codm, rcodm));
3705f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexTransformCreateDiscLabels(tr, rdm[r]));
3715f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetCoarseDM(rdm[r], cdm));
3725f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexSetRegularRefinement(rdm[r], PETSC_TRUE));
373012bc364SMatthew G. Knepley       if (rdm[r]) {
374012bc364SMatthew G. Knepley         ((DM_Plex *) (rdm[r])->data)->printFEM = ((DM_Plex *) dm->data)->printFEM;
375012bc364SMatthew G. Knepley         ((DM_Plex *) (rdm[r])->data)->printL2  = ((DM_Plex *) dm->data)->printL2;
3766d7c9049SMatthew G. Knepley       }
377012bc364SMatthew G. Knepley       cdm  = rdm[r];
3785f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexTransformDestroy(&tr));
3790d1cd5e0SMatthew G. Knepley     }
3800d1cd5e0SMatthew G. Knepley   } else {
3810d1cd5e0SMatthew G. Knepley     for (r = 0; r < nlevels; ++r) {
3825f80ce2aSJacob Faibussowitsch       CHKERRQ(DMRefine(cdm, PetscObjectComm((PetscObject) dm), &rdm[r]));
3835f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCopyDisc(cdm, rdm[r]));
3845f80ce2aSJacob Faibussowitsch       if (localized) CHKERRQ(DMLocalizeCoordinates(rdm[r]));
3855f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetCoarseDM(rdm[r], cdm));
386012bc364SMatthew G. Knepley       if (rdm[r]) {
387012bc364SMatthew G. Knepley         ((DM_Plex *) (rdm[r])->data)->printFEM = ((DM_Plex *) dm->data)->printFEM;
388012bc364SMatthew G. Knepley         ((DM_Plex *) (rdm[r])->data)->printL2  = ((DM_Plex *) dm->data)->printL2;
3896d7c9049SMatthew G. Knepley       }
390012bc364SMatthew G. Knepley       cdm  = rdm[r];
3910d1cd5e0SMatthew G. Knepley     }
3920d1cd5e0SMatthew G. Knepley   }
3930d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
3940d1cd5e0SMatthew G. Knepley }
395