xref: /petsc/src/dm/impls/plex/plexrefine.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 
22db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `DMPlexCreateTwoSidedProcessSF()`
23963fc26aSMatthew G. Knepley @*/
24*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateProcessSF(DM dm, PetscSF sfPoint, IS *processRanks, PetscSF *sfProcess)
25*d71ae5a4SJacob Faibussowitsch {
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);
37ad540459SPierre Jolivet   if (processRanks) PetscValidPointer(processRanks, 3);
38ad540459SPierre Jolivet   if (sfProcess) PetscValidPointer(sfProcess, 4);
399566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
409566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints));
419566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numLeaves, &ranks));
42ad540459SPierre Jolivet   for (l = 0; l < numLeaves; ++l) ranks[l] = remotePoints[l].rank;
439566063dSJacob Faibussowitsch   PetscCall(PetscSortRemoveDupsInt(&numLeaves, ranks));
449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numLeaves, &ranksNew));
459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numLeaves, &localPointsNew));
469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numLeaves, &remotePointsNew));
4775d3a19aSMatthew G. Knepley   for (l = 0; l < numLeaves; ++l) {
4875d3a19aSMatthew G. Knepley     ranksNew[l]              = ranks[l];
4975d3a19aSMatthew G. Knepley     localPointsNew[l]        = l;
5075d3a19aSMatthew G. Knepley     remotePointsNew[l].index = 0;
5175d3a19aSMatthew G. Knepley     remotePointsNew[l].rank  = ranksNew[l];
5275d3a19aSMatthew G. Knepley   }
539566063dSJacob Faibussowitsch   PetscCall(PetscFree(ranks));
549566063dSJacob Faibussowitsch   if (processRanks) PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), numLeaves, ranksNew, PETSC_OWN_POINTER, processRanks));
559566063dSJacob Faibussowitsch   else PetscCall(PetscFree(ranksNew));
56963fc26aSMatthew G. Knepley   if (sfProcess) {
579566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess));
589566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)*sfProcess, "Process SF"));
599566063dSJacob Faibussowitsch     PetscCall(PetscSFSetFromOptions(*sfProcess));
609566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*sfProcess, size, numLeaves, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER));
61963fc26aSMatthew G. Knepley   }
6275d3a19aSMatthew G. Knepley   PetscFunctionReturn(0);
6375d3a19aSMatthew G. Knepley }
6475d3a19aSMatthew G. Knepley 
652389894bSMatthew G. Knepley /*@
662389894bSMatthew G. Knepley   DMPlexCreateCoarsePointIS - Creates an IS covering the coarse DM chart with the fine points as data
672389894bSMatthew G. Knepley 
68012bc364SMatthew G. Knepley   Collective on dm
69012bc364SMatthew G. Knepley 
702389894bSMatthew G. Knepley   Input Parameter:
712389894bSMatthew G. Knepley . dm - The coarse DM
722389894bSMatthew G. Knepley 
732389894bSMatthew G. Knepley   Output Parameter:
742389894bSMatthew G. Knepley . fpointIS - The IS of all the fine points which exist in the original coarse mesh
752389894bSMatthew G. Knepley 
762389894bSMatthew G. Knepley   Level: developer
772389894bSMatthew G. Knepley 
78db781477SPatrick Sanan .seealso: `DMRefine()`, `DMPlexSetRefinementUniform()`, `DMPlexGetSubpointIS()`
792389894bSMatthew G. Knepley @*/
80*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateCoarsePointIS(DM dm, IS *fpointIS)
81*d71ae5a4SJacob Faibussowitsch {
82012bc364SMatthew G. Knepley   DMPlexTransform tr;
83412e9a14SMatthew G. Knepley   PetscInt       *fpoints;
84327c2912SStefano Zampini   PetscInt        pStart, pEnd, p, vStart, vEnd, v;
852389894bSMatthew G. Knepley 
862389894bSMatthew G. Knepley   PetscFunctionBegin;
879566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
889566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
899566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
909566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetUp(tr));
919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pEnd - pStart, &fpoints));
922389894bSMatthew G. Knepley   for (p = 0; p < pEnd - pStart; ++p) fpoints[p] = -1;
93412e9a14SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
94012bc364SMatthew G. Knepley     PetscInt vNew = -1; /* quiet overzealous may be used uninitialized check */
95327c2912SStefano Zampini 
969566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformGetTargetPoint(tr, DM_POLYTOPE_POINT, DM_POLYTOPE_POINT, p, 0, &vNew));
97412e9a14SMatthew G. Knepley     fpoints[v - pStart] = vNew;
982389894bSMatthew G. Knepley   }
999566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformDestroy(&tr));
1009566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, fpoints, PETSC_OWN_POINTER, fpointIS));
1012389894bSMatthew G. Knepley   PetscFunctionReturn(0);
1022389894bSMatthew G. Knepley }
1032389894bSMatthew G. Knepley 
104012bc364SMatthew G. Knepley /*@C
105012bc364SMatthew G. Knepley   DMPlexSetTransformType - Set the transform type for uniform refinement
106012bc364SMatthew G. Knepley 
107012bc364SMatthew G. Knepley   Input Parameters:
108012bc364SMatthew G. Knepley + dm - The DM
109012bc364SMatthew G. Knepley - type - The transform type for uniform refinement
110012bc364SMatthew G. Knepley 
111012bc364SMatthew G. Knepley   Level: developer
112012bc364SMatthew G. Knepley 
113db781477SPatrick Sanan .seealso: `DMPlexTransformType`, `DMRefine()`, `DMPlexGetTransformType()`, `DMPlexSetRefinementUniform()`
114012bc364SMatthew G. Knepley @*/
115*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetTransformType(DM dm, DMPlexTransformType type)
116*d71ae5a4SJacob Faibussowitsch {
117012bc364SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
118012bc364SMatthew G. Knepley 
119012bc364SMatthew G. Knepley   PetscFunctionBegin;
120012bc364SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
121012bc364SMatthew G. Knepley   if (type) PetscValidCharPointer(type, 2);
1229566063dSJacob Faibussowitsch   PetscCall(PetscFree(mesh->transformType));
1239566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(type, &mesh->transformType));
124012bc364SMatthew G. Knepley   PetscFunctionReturn(0);
125012bc364SMatthew G. Knepley }
126012bc364SMatthew G. Knepley 
127012bc364SMatthew G. Knepley /*@C
128012bc364SMatthew G. Knepley   DMPlexGetTransformType - Retrieve the transform type for uniform refinement
129012bc364SMatthew G. Knepley 
130012bc364SMatthew G. Knepley   Input Parameter:
131012bc364SMatthew G. Knepley . dm - The DM
132012bc364SMatthew G. Knepley 
133012bc364SMatthew G. Knepley   Output Parameter:
134012bc364SMatthew G. Knepley . type - The transform type for uniform refinement
135012bc364SMatthew G. Knepley 
136012bc364SMatthew G. Knepley   Level: developer
137012bc364SMatthew G. Knepley 
138db781477SPatrick Sanan .seealso: `DMPlexTransformType`, `DMRefine()`, `DMPlexSetTransformType()`, `DMPlexGetRefinementUniform()`
139012bc364SMatthew G. Knepley @*/
140*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetTransformType(DM dm, DMPlexTransformType *type)
141*d71ae5a4SJacob Faibussowitsch {
142012bc364SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
143012bc364SMatthew G. Knepley 
144012bc364SMatthew G. Knepley   PetscFunctionBegin;
145012bc364SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
146012bc364SMatthew G. Knepley   PetscValidPointer(type, 2);
147012bc364SMatthew G. Knepley   *type = mesh->transformType;
148012bc364SMatthew G. Knepley   PetscFunctionReturn(0);
149012bc364SMatthew G. Knepley }
150012bc364SMatthew G. Knepley 
1510e2b6761SMatthew G. Knepley /*@
1520e2b6761SMatthew G. Knepley   DMPlexSetRefinementUniform - Set the flag for uniform refinement
1530e2b6761SMatthew G. Knepley 
1540e2b6761SMatthew G. Knepley   Input Parameters:
1550e2b6761SMatthew G. Knepley + dm - The DM
1560e2b6761SMatthew G. Knepley - refinementUniform - The flag for uniform refinement
1570e2b6761SMatthew G. Knepley 
1580e2b6761SMatthew G. Knepley   Level: developer
1590e2b6761SMatthew G. Knepley 
160db781477SPatrick Sanan .seealso: `DMRefine()`, `DMPlexGetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
1610e2b6761SMatthew G. Knepley @*/
162*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetRefinementUniform(DM dm, PetscBool refinementUniform)
163*d71ae5a4SJacob Faibussowitsch {
16475d3a19aSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
16575d3a19aSMatthew G. Knepley 
16675d3a19aSMatthew G. Knepley   PetscFunctionBegin;
167012bc364SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
16875d3a19aSMatthew G. Knepley   mesh->refinementUniform = refinementUniform;
16975d3a19aSMatthew G. Knepley   PetscFunctionReturn(0);
17075d3a19aSMatthew G. Knepley }
17175d3a19aSMatthew G. Knepley 
1720e2b6761SMatthew G. Knepley /*@
1730e2b6761SMatthew G. Knepley   DMPlexGetRefinementUniform - Retrieve the flag for uniform refinement
1740e2b6761SMatthew G. Knepley 
1750e2b6761SMatthew G. Knepley   Input Parameter:
1760e2b6761SMatthew G. Knepley . dm - The DM
1770e2b6761SMatthew G. Knepley 
1780e2b6761SMatthew G. Knepley   Output Parameter:
1790e2b6761SMatthew G. Knepley . refinementUniform - The flag for uniform refinement
1800e2b6761SMatthew G. Knepley 
1810e2b6761SMatthew G. Knepley   Level: developer
1820e2b6761SMatthew G. Knepley 
183db781477SPatrick Sanan .seealso: `DMRefine()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
1840e2b6761SMatthew G. Knepley @*/
185*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetRefinementUniform(DM dm, PetscBool *refinementUniform)
186*d71ae5a4SJacob Faibussowitsch {
18775d3a19aSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
18875d3a19aSMatthew G. Knepley 
18975d3a19aSMatthew G. Knepley   PetscFunctionBegin;
190012bc364SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
191dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(refinementUniform, 2);
19275d3a19aSMatthew G. Knepley   *refinementUniform = mesh->refinementUniform;
19375d3a19aSMatthew G. Knepley   PetscFunctionReturn(0);
19475d3a19aSMatthew G. Knepley }
19575d3a19aSMatthew G. Knepley 
1960e2b6761SMatthew G. Knepley /*@
1970e2b6761SMatthew G. Knepley   DMPlexSetRefinementLimit - Set the maximum cell volume for refinement
1980e2b6761SMatthew G. Knepley 
1990e2b6761SMatthew G. Knepley   Input Parameters:
2000e2b6761SMatthew G. Knepley + dm - The DM
2010e2b6761SMatthew G. Knepley - refinementLimit - The maximum cell volume in the refined mesh
2020e2b6761SMatthew G. Knepley 
2030e2b6761SMatthew G. Knepley   Level: developer
2040e2b6761SMatthew G. Knepley 
205db781477SPatrick Sanan .seealso: `DMRefine()`, `DMPlexGetRefinementLimit()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`
2060e2b6761SMatthew G. Knepley @*/
207*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetRefinementLimit(DM dm, PetscReal refinementLimit)
208*d71ae5a4SJacob Faibussowitsch {
20975d3a19aSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
21075d3a19aSMatthew G. Knepley 
21175d3a19aSMatthew G. Knepley   PetscFunctionBegin;
21275d3a19aSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
21375d3a19aSMatthew G. Knepley   mesh->refinementLimit = refinementLimit;
21475d3a19aSMatthew G. Knepley   PetscFunctionReturn(0);
21575d3a19aSMatthew G. Knepley }
21675d3a19aSMatthew G. Knepley 
2170e2b6761SMatthew G. Knepley /*@
2180e2b6761SMatthew G. Knepley   DMPlexGetRefinementLimit - Retrieve the maximum cell volume for refinement
2190e2b6761SMatthew G. Knepley 
2200e2b6761SMatthew G. Knepley   Input Parameter:
2210e2b6761SMatthew G. Knepley . dm - The DM
2220e2b6761SMatthew G. Knepley 
2230e2b6761SMatthew G. Knepley   Output Parameter:
2240e2b6761SMatthew G. Knepley . refinementLimit - The maximum cell volume in the refined mesh
2250e2b6761SMatthew G. Knepley 
2260e2b6761SMatthew G. Knepley   Level: developer
2270e2b6761SMatthew G. Knepley 
228db781477SPatrick Sanan .seealso: `DMRefine()`, `DMPlexSetRefinementLimit()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`
2290e2b6761SMatthew G. Knepley @*/
230*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetRefinementLimit(DM dm, PetscReal *refinementLimit)
231*d71ae5a4SJacob Faibussowitsch {
23275d3a19aSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
23375d3a19aSMatthew G. Knepley 
23475d3a19aSMatthew G. Knepley   PetscFunctionBegin;
23575d3a19aSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
236dadcf809SJacob Faibussowitsch   PetscValidRealPointer(refinementLimit, 2);
23775d3a19aSMatthew G. Knepley   /* if (mesh->refinementLimit < 0) = getMaxVolume()/2.0; */
23875d3a19aSMatthew G. Knepley   *refinementLimit = mesh->refinementLimit;
23975d3a19aSMatthew G. Knepley   PetscFunctionReturn(0);
24075d3a19aSMatthew G. Knepley }
24175d3a19aSMatthew G. Knepley 
242b28003e6SMatthew G. Knepley /*@
243b28003e6SMatthew G. Knepley   DMPlexSetRefinementFunction - Set the function giving the maximum cell volume for refinement
244b28003e6SMatthew G. Knepley 
245b28003e6SMatthew G. Knepley   Input Parameters:
246b28003e6SMatthew G. Knepley + dm - The DM
247b28003e6SMatthew G. Knepley - refinementFunc - Function giving the maximum cell volume in the refined mesh
248b28003e6SMatthew G. Knepley 
249b28003e6SMatthew G. Knepley   Note: The calling sequence is refinementFunc(coords, limit)
250b28003e6SMatthew G. Knepley $ coords - Coordinates of the current point, usually a cell centroid
251b28003e6SMatthew G. Knepley $ limit  - The maximum cell volume for a cell containing this point
252b28003e6SMatthew G. Knepley 
253b28003e6SMatthew G. Knepley   Level: developer
254b28003e6SMatthew G. Knepley 
255db781477SPatrick Sanan .seealso: `DMRefine()`, `DMPlexGetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
256b28003e6SMatthew G. Knepley @*/
257*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetRefinementFunction(DM dm, PetscErrorCode (*refinementFunc)(const PetscReal[], PetscReal *))
258*d71ae5a4SJacob Faibussowitsch {
259b28003e6SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
260b28003e6SMatthew G. Knepley 
261b28003e6SMatthew G. Knepley   PetscFunctionBegin;
262b28003e6SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
263b28003e6SMatthew G. Knepley   mesh->refinementFunc = refinementFunc;
264b28003e6SMatthew G. Knepley   PetscFunctionReturn(0);
265b28003e6SMatthew G. Knepley }
266b28003e6SMatthew G. Knepley 
267b28003e6SMatthew G. Knepley /*@
268b28003e6SMatthew G. Knepley   DMPlexGetRefinementFunction - Get the function giving the maximum cell volume for refinement
269b28003e6SMatthew G. Knepley 
270b28003e6SMatthew G. Knepley   Input Parameter:
271b28003e6SMatthew G. Knepley . dm - The DM
272b28003e6SMatthew G. Knepley 
273b28003e6SMatthew G. Knepley   Output Parameter:
274b28003e6SMatthew G. Knepley . refinementFunc - Function giving the maximum cell volume in the refined mesh
275b28003e6SMatthew G. Knepley 
276b28003e6SMatthew G. Knepley   Note: The calling sequence is refinementFunc(coords, limit)
277b28003e6SMatthew G. Knepley $ coords - Coordinates of the current point, usually a cell centroid
278b28003e6SMatthew G. Knepley $ limit  - The maximum cell volume for a cell containing this point
279b28003e6SMatthew G. Knepley 
280b28003e6SMatthew G. Knepley   Level: developer
281b28003e6SMatthew G. Knepley 
282db781477SPatrick Sanan .seealso: `DMRefine()`, `DMPlexSetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
283b28003e6SMatthew G. Knepley @*/
284*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetRefinementFunction(DM dm, PetscErrorCode (**refinementFunc)(const PetscReal[], PetscReal *))
285*d71ae5a4SJacob Faibussowitsch {
286b28003e6SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
287b28003e6SMatthew G. Knepley 
288b28003e6SMatthew G. Knepley   PetscFunctionBegin;
289b28003e6SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
290b28003e6SMatthew G. Knepley   PetscValidPointer(refinementFunc, 2);
291b28003e6SMatthew G. Knepley   *refinementFunc = mesh->refinementFunc;
292b28003e6SMatthew G. Knepley   PetscFunctionReturn(0);
293b28003e6SMatthew G. Knepley }
294b28003e6SMatthew G. Knepley 
295*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *rdm)
296*d71ae5a4SJacob Faibussowitsch {
297492b8470SStefano Zampini   PetscBool isUniform;
2980d1cd5e0SMatthew G. Knepley 
2990d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
3009566063dSJacob Faibussowitsch   PetscCall(DMPlexGetRefinementUniform(dm, &isUniform));
3019566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-initref_dm_view"));
3020d1cd5e0SMatthew G. Knepley   if (isUniform) {
303012bc364SMatthew G. Knepley     DMPlexTransform     tr;
30451a74b61SMatthew G. Knepley     DM                  cdm, rcdm;
305012bc364SMatthew G. Knepley     DMPlexTransformType trType;
306012bc364SMatthew G. Knepley     const char         *prefix;
307012bc364SMatthew G. Knepley     PetscOptions        options;
3080d1cd5e0SMatthew G. Knepley 
3099566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
3109566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformSetDM(tr, dm));
3119566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransformType(dm, &trType));
3129566063dSJacob Faibussowitsch     if (trType) PetscCall(DMPlexTransformSetType(tr, trType));
3139566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
3149566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
3159566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptions((PetscObject)dm, &options));
3169566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptions((PetscObject)tr, options));
3179566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformSetFromOptions(tr));
3189566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL));
3199566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformSetUp(tr));
3209566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_transform_view"));
3219566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformApply(tr, dm, rdm));
3229566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRegularRefinement(*rdm, PETSC_TRUE));
3239566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, *rdm));
3249566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm, &cdm));
3259566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(*rdm, &rcdm));
3269566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(cdm, rcdm));
3279566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformCreateDiscLabels(tr, *rdm));
3289566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformDestroy(&tr));
3290d1cd5e0SMatthew G. Knepley   } else {
3309566063dSJacob Faibussowitsch     PetscCall(DMPlexRefine_Internal(dm, NULL, NULL, NULL, rdm));
3310d1cd5e0SMatthew G. Knepley   }
332012bc364SMatthew G. Knepley   if (*rdm) {
333012bc364SMatthew G. Knepley     ((DM_Plex *)(*rdm)->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
334012bc364SMatthew G. Knepley     ((DM_Plex *)(*rdm)->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
3356d7c9049SMatthew G. Knepley   }
3369566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-postref_dm_view"));
3370d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
3380d1cd5e0SMatthew G. Knepley }
3390d1cd5e0SMatthew G. Knepley 
340*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM rdm[])
341*d71ae5a4SJacob Faibussowitsch {
3420d1cd5e0SMatthew G. Knepley   DM        cdm = dm;
3430d1cd5e0SMatthew G. Knepley   PetscInt  r;
3440d1cd5e0SMatthew G. Knepley   PetscBool isUniform, localized;
3450d1cd5e0SMatthew G. Knepley 
3460d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
3479566063dSJacob Faibussowitsch   PetscCall(DMPlexGetRefinementUniform(dm, &isUniform));
3489566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalized(dm, &localized));
3490d1cd5e0SMatthew G. Knepley   if (isUniform) {
3500d1cd5e0SMatthew G. Knepley     for (r = 0; r < nlevels; ++r) {
351012bc364SMatthew G. Knepley       DMPlexTransform tr;
35251a74b61SMatthew G. Knepley       DM              codm, rcodm;
353012bc364SMatthew G. Knepley       const char     *prefix;
3540d1cd5e0SMatthew G. Knepley 
3559566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)cdm), &tr));
3569566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdm, &prefix));
3579566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
3589566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformSetDM(tr, cdm));
3599566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformSetFromOptions(tr));
3609566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformSetUp(tr));
3619566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformApply(tr, cdm, &rdm[r]));
3629566063dSJacob Faibussowitsch       PetscCall(DMSetCoarsenLevel(rdm[r], cdm->leveldown));
3639566063dSJacob Faibussowitsch       PetscCall(DMSetRefineLevel(rdm[r], cdm->levelup + 1));
3649566063dSJacob Faibussowitsch       PetscCall(DMCopyDisc(cdm, rdm[r]));
3659566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinateDM(dm, &codm));
3669566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinateDM(rdm[r], &rcodm));
3679566063dSJacob Faibussowitsch       PetscCall(DMCopyDisc(codm, rcodm));
3689566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformCreateDiscLabels(tr, rdm[r]));
3699566063dSJacob Faibussowitsch       PetscCall(DMSetCoarseDM(rdm[r], cdm));
3709566063dSJacob Faibussowitsch       PetscCall(DMPlexSetRegularRefinement(rdm[r], PETSC_TRUE));
371012bc364SMatthew G. Knepley       if (rdm[r]) {
372012bc364SMatthew G. Knepley         ((DM_Plex *)(rdm[r])->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
373012bc364SMatthew G. Knepley         ((DM_Plex *)(rdm[r])->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
3746d7c9049SMatthew G. Knepley       }
375012bc364SMatthew G. Knepley       cdm = rdm[r];
3769566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformDestroy(&tr));
3770d1cd5e0SMatthew G. Knepley     }
3780d1cd5e0SMatthew G. Knepley   } else {
3790d1cd5e0SMatthew G. Knepley     for (r = 0; r < nlevels; ++r) {
3809566063dSJacob Faibussowitsch       PetscCall(DMRefine(cdm, PetscObjectComm((PetscObject)dm), &rdm[r]));
3819566063dSJacob Faibussowitsch       PetscCall(DMCopyDisc(cdm, rdm[r]));
3829566063dSJacob Faibussowitsch       if (localized) PetscCall(DMLocalizeCoordinates(rdm[r]));
3839566063dSJacob Faibussowitsch       PetscCall(DMSetCoarseDM(rdm[r], cdm));
384012bc364SMatthew G. Knepley       if (rdm[r]) {
385012bc364SMatthew G. Knepley         ((DM_Plex *)(rdm[r])->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
386012bc364SMatthew G. Knepley         ((DM_Plex *)(rdm[r])->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
3876d7c9049SMatthew G. Knepley       }
388012bc364SMatthew G. Knepley       cdm = rdm[r];
3890d1cd5e0SMatthew G. Knepley     }
3900d1cd5e0SMatthew G. Knepley   }
3910d1cd5e0SMatthew G. Knepley   PetscFunctionReturn(0);
3920d1cd5e0SMatthew G. Knepley }
393