xref: /petsc/src/dm/impls/plex/plexrefine.c (revision 3f5a071cd8576c6a48d5fe37e5d896863c904fc8)
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 
4cc4c1da9SBarry Smith #include <petscdmplextransform.h> /*I      "petscdmplextransform.h"   I*/
575d3a19aSMatthew G. Knepley #include <petscsf.h>
675d3a19aSMatthew G. Knepley 
7963fc26aSMatthew G. Knepley /*@
8a1cb98faSBarry Smith   DMPlexCreateProcessSF - Create an `PetscSF` which just has process connectivity
9963fc26aSMatthew G. Knepley 
1020f4b53cSBarry Smith   Collective
11963fc26aSMatthew G. Knepley 
12963fc26aSMatthew G. Knepley   Input Parameters:
13a1cb98faSBarry Smith + dm      - The `DM`
14a1cb98faSBarry Smith - sfPoint - The `PetscSF` which encodes point connectivity
15963fc26aSMatthew G. Knepley 
16963fc26aSMatthew G. Knepley   Output Parameters:
1720f4b53cSBarry Smith + processRanks - A list of process neighbors, or `NULL`
1820f4b53cSBarry Smith - sfProcess    - An `PetscSF` encoding the process connectivity, or `NULL`
19963fc26aSMatthew G. Knepley 
20963fc26aSMatthew G. Knepley   Level: developer
21963fc26aSMatthew G. Knepley 
221cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `PetscSF`, `PetscSFCreate()`, `DMPlexCreateTwoSidedProcessSF()`
23963fc26aSMatthew G. Knepley @*/
DMPlexCreateProcessSF(DM dm,PetscSF sfPoint,IS * processRanks,PetscSF * sfProcess)24d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateProcessSF(DM dm, PetscSF sfPoint, IS *processRanks, PetscSF *sfProcess)
25d71ae5a4SJacob 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;
316497c311SBarry Smith   PetscMPIInt       *ranks;
326497c311SBarry Smith   PetscInt          *ranksNew;
339852e123SBarry Smith   PetscMPIInt        size;
3475d3a19aSMatthew G. Knepley 
3575d3a19aSMatthew G. Knepley   PetscFunctionBegin;
36963fc26aSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
37963fc26aSMatthew G. Knepley   PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2);
384f572ea9SToby Isaac   if (processRanks) PetscAssertPointer(processRanks, 3);
394f572ea9SToby Isaac   if (sfProcess) PetscAssertPointer(sfProcess, 4);
409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
419566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints));
429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numLeaves, &ranks));
436497c311SBarry Smith   for (l = 0; l < numLeaves; ++l) ranks[l] = (PetscMPIInt)remotePoints[l].rank;
446497c311SBarry Smith   PetscCall(PetscSortRemoveDupsMPIInt(&numLeaves, ranks));
459566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numLeaves, &ranksNew));
469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numLeaves, &localPointsNew));
479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numLeaves, &remotePointsNew));
4875d3a19aSMatthew G. Knepley   for (l = 0; l < numLeaves; ++l) {
4975d3a19aSMatthew G. Knepley     ranksNew[l]              = ranks[l];
5075d3a19aSMatthew G. Knepley     localPointsNew[l]        = l;
5175d3a19aSMatthew G. Knepley     remotePointsNew[l].index = 0;
52835f2295SStefano Zampini     remotePointsNew[l].rank  = ranksNew[l];
5375d3a19aSMatthew G. Knepley   }
549566063dSJacob Faibussowitsch   PetscCall(PetscFree(ranks));
559566063dSJacob Faibussowitsch   if (processRanks) PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), numLeaves, ranksNew, PETSC_OWN_POINTER, processRanks));
569566063dSJacob Faibussowitsch   else PetscCall(PetscFree(ranksNew));
57963fc26aSMatthew G. Knepley   if (sfProcess) {
589566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess));
599566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)*sfProcess, "Process SF"));
609566063dSJacob Faibussowitsch     PetscCall(PetscSFSetFromOptions(*sfProcess));
619566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*sfProcess, size, numLeaves, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER));
62963fc26aSMatthew G. Knepley   }
633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6475d3a19aSMatthew G. Knepley }
6575d3a19aSMatthew G. Knepley 
662389894bSMatthew G. Knepley /*@
67a1cb98faSBarry Smith   DMPlexCreateCoarsePointIS - Creates an `IS` covering the coarse `DM` chart with the fine points as data
682389894bSMatthew G. Knepley 
6920f4b53cSBarry Smith   Collective
70012bc364SMatthew G. Knepley 
712389894bSMatthew G. Knepley   Input Parameter:
72a1cb98faSBarry Smith . dm - The coarse `DM`
732389894bSMatthew G. Knepley 
742389894bSMatthew G. Knepley   Output Parameter:
75a1cb98faSBarry Smith . fpointIS - The `IS` of all the fine points which exist in the original coarse mesh
762389894bSMatthew G. Knepley 
772389894bSMatthew G. Knepley   Level: developer
782389894bSMatthew G. Knepley 
791cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `IS`, `DMRefine()`, `DMPlexSetRefinementUniform()`, `DMPlexGetSubpointIS()`
802389894bSMatthew G. Knepley @*/
DMPlexCreateCoarsePointIS(DM dm,IS * fpointIS)81d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateCoarsePointIS(DM dm, IS *fpointIS)
82d71ae5a4SJacob Faibussowitsch {
83012bc364SMatthew G. Knepley   DMPlexTransform tr;
84412e9a14SMatthew G. Knepley   PetscInt       *fpoints;
85327c2912SStefano Zampini   PetscInt        pStart, pEnd, p, vStart, vEnd, v;
862389894bSMatthew G. Knepley 
872389894bSMatthew G. Knepley   PetscFunctionBegin;
889566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
899566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
909566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
919566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformSetUp(tr));
929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pEnd - pStart, &fpoints));
932389894bSMatthew G. Knepley   for (p = 0; p < pEnd - pStart; ++p) fpoints[p] = -1;
94412e9a14SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
95012bc364SMatthew G. Knepley     PetscInt vNew = -1; /* quiet overzealous may be used uninitialized check */
96327c2912SStefano Zampini 
979566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformGetTargetPoint(tr, DM_POLYTOPE_POINT, DM_POLYTOPE_POINT, p, 0, &vNew));
98412e9a14SMatthew G. Knepley     fpoints[v - pStart] = vNew;
992389894bSMatthew G. Knepley   }
1009566063dSJacob Faibussowitsch   PetscCall(DMPlexTransformDestroy(&tr));
1019566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, fpoints, PETSC_OWN_POINTER, fpointIS));
1023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1032389894bSMatthew G. Knepley }
1042389894bSMatthew G. Knepley 
105ce78bad3SBarry Smith /*@
106012bc364SMatthew G. Knepley   DMPlexSetTransformType - Set the transform type for uniform refinement
107012bc364SMatthew G. Knepley 
108012bc364SMatthew G. Knepley   Input Parameters:
109a1cb98faSBarry Smith + dm   - The `DM`
110012bc364SMatthew G. Knepley - type - The transform type for uniform refinement
111012bc364SMatthew G. Knepley 
112012bc364SMatthew G. Knepley   Level: developer
113012bc364SMatthew G. Knepley 
1141cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransformType`, `DMRefine()`, `DMPlexGetTransformType()`, `DMPlexSetRefinementUniform()`
115012bc364SMatthew G. Knepley @*/
DMPlexSetTransformType(DM dm,DMPlexTransformType type)116d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetTransformType(DM dm, DMPlexTransformType type)
117d71ae5a4SJacob Faibussowitsch {
118012bc364SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
119012bc364SMatthew G. Knepley 
120012bc364SMatthew G. Knepley   PetscFunctionBegin;
121012bc364SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
1224f572ea9SToby Isaac   if (type) PetscAssertPointer(type, 2);
1239566063dSJacob Faibussowitsch   PetscCall(PetscFree(mesh->transformType));
1249566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(type, &mesh->transformType));
1253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
126012bc364SMatthew G. Knepley }
127012bc364SMatthew G. Knepley 
128012bc364SMatthew G. Knepley /*@C
129012bc364SMatthew G. Knepley   DMPlexGetTransformType - Retrieve the transform type for uniform refinement
130012bc364SMatthew G. Knepley 
131012bc364SMatthew G. Knepley   Input Parameter:
132a1cb98faSBarry Smith . dm - The `DM`
133012bc364SMatthew G. Knepley 
134012bc364SMatthew G. Knepley   Output Parameter:
135012bc364SMatthew G. Knepley . type - The transform type for uniform refinement
136012bc364SMatthew G. Knepley 
137012bc364SMatthew G. Knepley   Level: developer
138012bc364SMatthew G. Knepley 
1391cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexTransformType`, `DMRefine()`, `DMPlexSetTransformType()`, `DMPlexGetRefinementUniform()`
140012bc364SMatthew G. Knepley @*/
DMPlexGetTransformType(DM dm,DMPlexTransformType * type)141d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetTransformType(DM dm, DMPlexTransformType *type)
142d71ae5a4SJacob Faibussowitsch {
143012bc364SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
144012bc364SMatthew G. Knepley 
145012bc364SMatthew G. Knepley   PetscFunctionBegin;
146012bc364SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
1474f572ea9SToby Isaac   PetscAssertPointer(type, 2);
148012bc364SMatthew G. Knepley   *type = mesh->transformType;
1493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
150012bc364SMatthew G. Knepley }
151012bc364SMatthew G. Knepley 
DMPlexSetTransform(DM dm,DMPlexTransform tr)15261f058f9SMatthew G. Knepley PetscErrorCode DMPlexSetTransform(DM dm, DMPlexTransform tr)
15361f058f9SMatthew G. Knepley {
15461f058f9SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
15561f058f9SMatthew G. Knepley 
15661f058f9SMatthew G. Knepley   PetscFunctionBegin;
15761f058f9SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
15861f058f9SMatthew G. Knepley   if (tr) PetscValidHeaderSpecific(tr, DMPLEXTRANSFORM_CLASSID, 2);
15961f058f9SMatthew G. Knepley   PetscCall(PetscObjectReference((PetscObject)tr));
16061f058f9SMatthew G. Knepley   PetscCall(DMPlexTransformDestroy(&mesh->transform));
16161f058f9SMatthew G. Knepley   mesh->transform = tr;
16261f058f9SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
16361f058f9SMatthew G. Knepley }
16461f058f9SMatthew G. Knepley 
DMPlexGetTransform(DM dm,DMPlexTransform * tr)16561f058f9SMatthew G. Knepley PetscErrorCode DMPlexGetTransform(DM dm, DMPlexTransform *tr)
16661f058f9SMatthew G. Knepley {
16761f058f9SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
16861f058f9SMatthew G. Knepley 
16961f058f9SMatthew G. Knepley   PetscFunctionBegin;
17061f058f9SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
17161f058f9SMatthew G. Knepley   PetscAssertPointer(tr, 2);
17261f058f9SMatthew G. Knepley   *tr = mesh->transform;
17361f058f9SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
17461f058f9SMatthew G. Knepley }
17561f058f9SMatthew G. Knepley 
DMPlexSetSaveTransform(DM dm,PetscBool save)17661f058f9SMatthew G. Knepley PetscErrorCode DMPlexSetSaveTransform(DM dm, PetscBool save)
17761f058f9SMatthew G. Knepley {
17861f058f9SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
17961f058f9SMatthew G. Knepley 
18061f058f9SMatthew G. Knepley   PetscFunctionBegin;
18161f058f9SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
18261f058f9SMatthew G. Knepley   mesh->saveTransform = save;
18361f058f9SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
18461f058f9SMatthew G. Knepley }
18561f058f9SMatthew G. Knepley 
DMPlexGetSaveTransform(DM dm,PetscBool * save)18661f058f9SMatthew G. Knepley PetscErrorCode DMPlexGetSaveTransform(DM dm, PetscBool *save)
18761f058f9SMatthew G. Knepley {
18861f058f9SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
18961f058f9SMatthew G. Knepley 
19061f058f9SMatthew G. Knepley   PetscFunctionBegin;
19161f058f9SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
19261f058f9SMatthew G. Knepley   PetscAssertPointer(save, 2);
19361f058f9SMatthew G. Knepley   *save = mesh->saveTransform;
19461f058f9SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
19561f058f9SMatthew G. Knepley }
19661f058f9SMatthew G. Knepley 
1970e2b6761SMatthew G. Knepley /*@
1980e2b6761SMatthew G. Knepley   DMPlexSetRefinementUniform - Set the flag for uniform refinement
1990e2b6761SMatthew G. Knepley 
2000e2b6761SMatthew G. Knepley   Input Parameters:
201a1cb98faSBarry Smith + dm                - The `DM`
2020e2b6761SMatthew G. Knepley - refinementUniform - The flag for uniform refinement
2030e2b6761SMatthew G. Knepley 
2040e2b6761SMatthew G. Knepley   Level: developer
2050e2b6761SMatthew G. Knepley 
2061cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
2070e2b6761SMatthew G. Knepley @*/
DMPlexSetRefinementUniform(DM dm,PetscBool refinementUniform)208d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetRefinementUniform(DM dm, PetscBool refinementUniform)
209d71ae5a4SJacob Faibussowitsch {
21075d3a19aSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
21175d3a19aSMatthew G. Knepley 
21275d3a19aSMatthew G. Knepley   PetscFunctionBegin;
213012bc364SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
21475d3a19aSMatthew G. Knepley   mesh->refinementUniform = refinementUniform;
2153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21675d3a19aSMatthew G. Knepley }
21775d3a19aSMatthew G. Knepley 
2180e2b6761SMatthew G. Knepley /*@
2190e2b6761SMatthew G. Knepley   DMPlexGetRefinementUniform - Retrieve the flag for uniform refinement
2200e2b6761SMatthew G. Knepley 
2210e2b6761SMatthew G. Knepley   Input Parameter:
222a1cb98faSBarry Smith . dm - The `DM`
2230e2b6761SMatthew G. Knepley 
2240e2b6761SMatthew G. Knepley   Output Parameter:
2250e2b6761SMatthew G. Knepley . refinementUniform - The flag for uniform refinement
2260e2b6761SMatthew G. Knepley 
2270e2b6761SMatthew G. Knepley   Level: developer
2280e2b6761SMatthew G. Knepley 
2291cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
2300e2b6761SMatthew G. Knepley @*/
DMPlexGetRefinementUniform(DM dm,PetscBool * refinementUniform)231d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetRefinementUniform(DM dm, PetscBool *refinementUniform)
232d71ae5a4SJacob Faibussowitsch {
23375d3a19aSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
23475d3a19aSMatthew G. Knepley 
23575d3a19aSMatthew G. Knepley   PetscFunctionBegin;
236012bc364SMatthew G. Knepley   PetscValidHeaderSpecificType(dm, DM_CLASSID, 1, DMPLEX);
2374f572ea9SToby Isaac   PetscAssertPointer(refinementUniform, 2);
23875d3a19aSMatthew G. Knepley   *refinementUniform = mesh->refinementUniform;
2393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
24075d3a19aSMatthew G. Knepley }
24175d3a19aSMatthew G. Knepley 
2420e2b6761SMatthew G. Knepley /*@
2430e2b6761SMatthew G. Knepley   DMPlexSetRefinementLimit - Set the maximum cell volume for refinement
2440e2b6761SMatthew G. Knepley 
2450e2b6761SMatthew G. Knepley   Input Parameters:
246a1cb98faSBarry Smith + dm              - The `DM`
2470e2b6761SMatthew G. Knepley - refinementLimit - The maximum cell volume in the refined mesh
2480e2b6761SMatthew G. Knepley 
2490e2b6761SMatthew G. Knepley   Level: developer
2500e2b6761SMatthew G. Knepley 
2511cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementLimit()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`
2520e2b6761SMatthew G. Knepley @*/
DMPlexSetRefinementLimit(DM dm,PetscReal refinementLimit)253d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSetRefinementLimit(DM dm, PetscReal refinementLimit)
254d71ae5a4SJacob Faibussowitsch {
25575d3a19aSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
25675d3a19aSMatthew G. Knepley 
25775d3a19aSMatthew G. Knepley   PetscFunctionBegin;
25875d3a19aSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
25975d3a19aSMatthew G. Knepley   mesh->refinementLimit = refinementLimit;
2603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26175d3a19aSMatthew G. Knepley }
26275d3a19aSMatthew G. Knepley 
2630e2b6761SMatthew G. Knepley /*@
2640e2b6761SMatthew G. Knepley   DMPlexGetRefinementLimit - Retrieve the maximum cell volume for refinement
2650e2b6761SMatthew G. Knepley 
2660e2b6761SMatthew G. Knepley   Input Parameter:
267a1cb98faSBarry Smith . dm - The `DM`
2680e2b6761SMatthew G. Knepley 
2690e2b6761SMatthew G. Knepley   Output Parameter:
2700e2b6761SMatthew G. Knepley . refinementLimit - The maximum cell volume in the refined mesh
2710e2b6761SMatthew G. Knepley 
2720e2b6761SMatthew G. Knepley   Level: developer
2730e2b6761SMatthew G. Knepley 
2741cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementLimit()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`
2750e2b6761SMatthew G. Knepley @*/
DMPlexGetRefinementLimit(DM dm,PetscReal * refinementLimit)276d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexGetRefinementLimit(DM dm, PetscReal *refinementLimit)
277d71ae5a4SJacob Faibussowitsch {
27875d3a19aSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
27975d3a19aSMatthew G. Knepley 
28075d3a19aSMatthew G. Knepley   PetscFunctionBegin;
28175d3a19aSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2824f572ea9SToby Isaac   PetscAssertPointer(refinementLimit, 2);
28375d3a19aSMatthew G. Knepley   /* if (mesh->refinementLimit < 0) = getMaxVolume()/2.0; */
28475d3a19aSMatthew G. Knepley   *refinementLimit = mesh->refinementLimit;
2853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28675d3a19aSMatthew G. Knepley }
28775d3a19aSMatthew G. Knepley 
28860225df5SJacob Faibussowitsch /*@C
289b28003e6SMatthew G. Knepley   DMPlexSetRefinementFunction - Set the function giving the maximum cell volume for refinement
290b28003e6SMatthew G. Knepley 
291b28003e6SMatthew G. Knepley   Input Parameters:
292a1cb98faSBarry Smith + dm             - The `DM`
293b28003e6SMatthew G. Knepley - refinementFunc - Function giving the maximum cell volume in the refined mesh
294b28003e6SMatthew G. Knepley 
29520f4b53cSBarry Smith   Calling Sequence of `refinementFunc`:
296a1cb98faSBarry Smith + coords - Coordinates of the current point, usually a cell centroid
297a1cb98faSBarry Smith - limit  - The maximum cell volume for a cell containing this point
298b28003e6SMatthew G. Knepley 
299b28003e6SMatthew G. Knepley   Level: developer
300b28003e6SMatthew G. Knepley 
3011cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexGetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
302b28003e6SMatthew G. Knepley @*/
DMPlexSetRefinementFunction(DM dm,PetscErrorCode (* refinementFunc)(const PetscReal coords[],PetscReal * limit))303a4e35b19SJacob Faibussowitsch PetscErrorCode DMPlexSetRefinementFunction(DM dm, PetscErrorCode (*refinementFunc)(const PetscReal coords[], PetscReal *limit))
304d71ae5a4SJacob Faibussowitsch {
305b28003e6SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
306b28003e6SMatthew G. Knepley 
307b28003e6SMatthew G. Knepley   PetscFunctionBegin;
308b28003e6SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
309b28003e6SMatthew G. Knepley   mesh->refinementFunc = refinementFunc;
3103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
311b28003e6SMatthew G. Knepley }
312b28003e6SMatthew G. Knepley 
31360225df5SJacob Faibussowitsch /*@C
314b28003e6SMatthew G. Knepley   DMPlexGetRefinementFunction - Get the function giving the maximum cell volume for refinement
315b28003e6SMatthew G. Knepley 
316b28003e6SMatthew G. Knepley   Input Parameter:
31720f4b53cSBarry Smith . dm - The `DM`
318b28003e6SMatthew G. Knepley 
319b28003e6SMatthew G. Knepley   Output Parameter:
320b28003e6SMatthew G. Knepley . refinementFunc - Function giving the maximum cell volume in the refined mesh
321b28003e6SMatthew G. Knepley 
32220f4b53cSBarry Smith   Calling Sequence of `refinementFunc`:
323a1cb98faSBarry Smith + coords - Coordinates of the current point, usually a cell centroid
324a1cb98faSBarry Smith - limit  - The maximum cell volume for a cell containing this point
325b28003e6SMatthew G. Knepley 
326b28003e6SMatthew G. Knepley   Level: developer
327b28003e6SMatthew G. Knepley 
3281cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMRefine()`, `DMPlexSetRefinementFunction()`, `DMPlexGetRefinementUniform()`, `DMPlexSetRefinementUniform()`, `DMPlexGetRefinementLimit()`, `DMPlexSetRefinementLimit()`
329b28003e6SMatthew G. Knepley @*/
DMPlexGetRefinementFunction(DM dm,PetscErrorCode (** refinementFunc)(const PetscReal coords[],PetscReal * limit))330a4e35b19SJacob Faibussowitsch PetscErrorCode DMPlexGetRefinementFunction(DM dm, PetscErrorCode (**refinementFunc)(const PetscReal coords[], PetscReal *limit))
331d71ae5a4SJacob Faibussowitsch {
332b28003e6SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *)dm->data;
333b28003e6SMatthew G. Knepley 
334b28003e6SMatthew G. Knepley   PetscFunctionBegin;
335b28003e6SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3364f572ea9SToby Isaac   PetscAssertPointer(refinementFunc, 2);
337b28003e6SMatthew G. Knepley   *refinementFunc = mesh->refinementFunc;
3383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
339b28003e6SMatthew G. Knepley }
340b28003e6SMatthew G. Knepley 
DMRefine_Plex(DM dm,MPI_Comm comm,DM * rdm)341d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefine_Plex(DM dm, MPI_Comm comm, DM *rdm)
342d71ae5a4SJacob Faibussowitsch {
343492b8470SStefano Zampini   PetscBool isUniform;
3440d1cd5e0SMatthew G. Knepley 
3450d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
3469566063dSJacob Faibussowitsch   PetscCall(DMPlexGetRefinementUniform(dm, &isUniform));
3479566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-initref_dm_view"));
3480d1cd5e0SMatthew G. Knepley   if (isUniform) {
349012bc364SMatthew G. Knepley     DMPlexTransform     tr;
35051a74b61SMatthew G. Knepley     DM                  cdm, rcdm;
351012bc364SMatthew G. Knepley     DMPlexTransformType trType;
352012bc364SMatthew G. Knepley     const char         *prefix;
353012bc364SMatthew G. Knepley     PetscOptions        options;
354e44f6aebSMatthew G. Knepley     PetscInt            cDegree;
3559171e0a4SMatthew G. Knepley     PetscBool           useCeed, save;
3560d1cd5e0SMatthew G. Knepley 
3579566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), &tr));
3589566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformSetDM(tr, dm));
3599566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransformType(dm, &trType));
3609566063dSJacob Faibussowitsch     if (trType) PetscCall(DMPlexTransformSetType(tr, trType));
3619566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
3629566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
3639566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptions((PetscObject)dm, &options));
3649566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptions((PetscObject)tr, options));
3659566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformSetFromOptions(tr));
3669566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptions((PetscObject)tr, NULL));
3679566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformSetUp(tr));
3689566063dSJacob Faibussowitsch     PetscCall(PetscObjectViewFromOptions((PetscObject)tr, NULL, "-dm_plex_transform_view"));
3699566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformApply(tr, dm, rdm));
3709566063dSJacob Faibussowitsch     PetscCall(DMPlexSetRegularRefinement(*rdm, PETSC_TRUE));
371d2b2dc1eSMatthew G. Knepley     PetscCall(DMPlexGetUseCeed(dm, &useCeed));
372d2b2dc1eSMatthew G. Knepley     PetscCall(DMPlexSetUseCeed(*rdm, useCeed));
373f4f49eeaSPierre Jolivet     PetscCall(DMSetMatType(*rdm, dm->mattype));
3749566063dSJacob Faibussowitsch     PetscCall(DMCopyDisc(dm, *rdm));
3759566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm, &cdm));
3769566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(*rdm, &rcdm));
377e44f6aebSMatthew G. Knepley     PetscCall(DMGetCoordinateDegree_Internal(dm, &cDegree));
37838cc584fSMatthew G. Knepley     {
37938cc584fSMatthew G. Knepley       PetscDS cds, rcds;
38038cc584fSMatthew G. Knepley 
381*e65c294aSksagiyam       PetscCall(DMPlexCreateCoordinateSpace(*rdm, cDegree, PETSC_FALSE, PETSC_TRUE));
382e44f6aebSMatthew G. Knepley       PetscCall(DMGetCoordinateDM(*rdm, &rcdm));
38338cc584fSMatthew G. Knepley       PetscCall(DMGetDS(cdm, &cds));
38438cc584fSMatthew G. Knepley       PetscCall(DMGetDS(rcdm, &rcds));
38538cc584fSMatthew G. Knepley       PetscCall(PetscDSCopyConstants(cds, rcds));
386e44f6aebSMatthew G. Knepley     }
387d2b2dc1eSMatthew G. Knepley     PetscCall(DMPlexGetUseCeed(cdm, &useCeed));
388d2b2dc1eSMatthew G. Knepley     PetscCall(DMPlexSetUseCeed(rcdm, useCeed));
389d2b2dc1eSMatthew G. Knepley     if (useCeed) {
390d2b2dc1eSMatthew G. Knepley       PetscCall(DMPlexSetUseMatClosurePermutation(rcdm, PETSC_FALSE));
391d2b2dc1eSMatthew G. Knepley       PetscCall(DMUseTensorOrder(rcdm, PETSC_TRUE));
392d2b2dc1eSMatthew G. Knepley     }
3939566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformCreateDiscLabels(tr, *rdm));
3949171e0a4SMatthew G. Knepley     PetscCall(DMPlexGetSaveTransform(dm, &save));
3959171e0a4SMatthew G. Knepley     if (save) PetscCall(DMPlexSetTransform(*rdm, tr));
3969566063dSJacob Faibussowitsch     PetscCall(DMPlexTransformDestroy(&tr));
3970d1cd5e0SMatthew G. Knepley   } else {
3989566063dSJacob Faibussowitsch     PetscCall(DMPlexRefine_Internal(dm, NULL, NULL, NULL, rdm));
3990d1cd5e0SMatthew G. Knepley   }
400012bc364SMatthew G. Knepley   if (*rdm) {
401012bc364SMatthew G. Knepley     ((DM_Plex *)(*rdm)->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
402012bc364SMatthew G. Knepley     ((DM_Plex *)(*rdm)->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
4036d7c9049SMatthew G. Knepley   }
404e01335dfSMatthew G. Knepley   PetscCall(DMViewFromOptions(*rdm, NULL, "-postref_dm_view"));
4053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4060d1cd5e0SMatthew G. Knepley }
4070d1cd5e0SMatthew G. Knepley 
DMRefineHierarchy_Plex(DM dm,PetscInt nlevels,DM rdm[])408d71ae5a4SJacob Faibussowitsch PetscErrorCode DMRefineHierarchy_Plex(DM dm, PetscInt nlevels, DM rdm[])
409d71ae5a4SJacob Faibussowitsch {
4100d1cd5e0SMatthew G. Knepley   DM        cdm = dm;
4110d1cd5e0SMatthew G. Knepley   PetscInt  r;
412d2b2dc1eSMatthew G. Knepley   PetscBool isUniform, localized, useCeed;
4130d1cd5e0SMatthew G. Knepley 
4140d1cd5e0SMatthew G. Knepley   PetscFunctionBegin;
4159566063dSJacob Faibussowitsch   PetscCall(DMPlexGetRefinementUniform(dm, &isUniform));
4169566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalized(dm, &localized));
4170d1cd5e0SMatthew G. Knepley   if (isUniform) {
4180d1cd5e0SMatthew G. Knepley     for (r = 0; r < nlevels; ++r) {
419012bc364SMatthew G. Knepley       DMPlexTransform tr;
42051a74b61SMatthew G. Knepley       DM              codm, rcodm;
421012bc364SMatthew G. Knepley       const char     *prefix;
4220d1cd5e0SMatthew G. Knepley 
4239566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)cdm), &tr));
4249566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetOptionsPrefix((PetscObject)cdm, &prefix));
4259566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)tr, prefix));
4269566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformSetDM(tr, cdm));
4279566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformSetFromOptions(tr));
4289566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformSetUp(tr));
4299566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformApply(tr, cdm, &rdm[r]));
4309566063dSJacob Faibussowitsch       PetscCall(DMSetCoarsenLevel(rdm[r], cdm->leveldown));
4319566063dSJacob Faibussowitsch       PetscCall(DMSetRefineLevel(rdm[r], cdm->levelup + 1));
432bb393918Smarkadams4       PetscCall(DMSetMatType(rdm[r], dm->mattype));
433d2b2dc1eSMatthew G. Knepley       PetscCall(DMPlexGetUseCeed(dm, &useCeed));
434d2b2dc1eSMatthew G. Knepley       PetscCall(DMPlexSetUseCeed(rdm[r], useCeed));
4359566063dSJacob Faibussowitsch       PetscCall(DMCopyDisc(cdm, rdm[r]));
4369566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinateDM(dm, &codm));
4379566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinateDM(rdm[r], &rcodm));
4389566063dSJacob Faibussowitsch       PetscCall(DMCopyDisc(codm, rcodm));
439d2b2dc1eSMatthew G. Knepley       PetscCall(DMPlexGetUseCeed(codm, &useCeed));
440d2b2dc1eSMatthew G. Knepley       PetscCall(DMPlexSetUseCeed(rcodm, useCeed));
441d2b2dc1eSMatthew G. Knepley       if (useCeed) {
442d2b2dc1eSMatthew G. Knepley         PetscCall(DMPlexSetUseMatClosurePermutation(rcodm, PETSC_FALSE));
443d2b2dc1eSMatthew G. Knepley         PetscCall(DMUseTensorOrder(rcodm, PETSC_TRUE));
444d2b2dc1eSMatthew G. Knepley       }
4459566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformCreateDiscLabels(tr, rdm[r]));
4469566063dSJacob Faibussowitsch       PetscCall(DMSetCoarseDM(rdm[r], cdm));
4479566063dSJacob Faibussowitsch       PetscCall(DMPlexSetRegularRefinement(rdm[r], PETSC_TRUE));
448012bc364SMatthew G. Knepley       if (rdm[r]) {
449f4f49eeaSPierre Jolivet         ((DM_Plex *)rdm[r]->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
450f4f49eeaSPierre Jolivet         ((DM_Plex *)rdm[r]->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
4516d7c9049SMatthew G. Knepley       }
452012bc364SMatthew G. Knepley       cdm = rdm[r];
4539566063dSJacob Faibussowitsch       PetscCall(DMPlexTransformDestroy(&tr));
4540d1cd5e0SMatthew G. Knepley     }
4550d1cd5e0SMatthew G. Knepley   } else {
4560d1cd5e0SMatthew G. Knepley     for (r = 0; r < nlevels; ++r) {
4579566063dSJacob Faibussowitsch       PetscCall(DMRefine(cdm, PetscObjectComm((PetscObject)dm), &rdm[r]));
458d2b2dc1eSMatthew G. Knepley       PetscCall(DMPlexGetUseCeed(dm, &useCeed));
459d2b2dc1eSMatthew G. Knepley       PetscCall(DMPlexSetUseCeed(rdm[r], useCeed));
4609566063dSJacob Faibussowitsch       PetscCall(DMCopyDisc(cdm, rdm[r]));
4619566063dSJacob Faibussowitsch       if (localized) PetscCall(DMLocalizeCoordinates(rdm[r]));
4629566063dSJacob Faibussowitsch       PetscCall(DMSetCoarseDM(rdm[r], cdm));
463012bc364SMatthew G. Knepley       if (rdm[r]) {
464f4f49eeaSPierre Jolivet         ((DM_Plex *)rdm[r]->data)->printFEM = ((DM_Plex *)dm->data)->printFEM;
465f4f49eeaSPierre Jolivet         ((DM_Plex *)rdm[r]->data)->printL2  = ((DM_Plex *)dm->data)->printL2;
4666d7c9049SMatthew G. Knepley       }
467012bc364SMatthew G. Knepley       cdm = rdm[r];
4680d1cd5e0SMatthew G. Knepley     }
4690d1cd5e0SMatthew G. Knepley   }
4703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4710d1cd5e0SMatthew G. Knepley }
472