xref: /petsc/src/dm/impls/plex/tutorials/ex10.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1adcf2cb4SMatthew G. Knepley static char help[] = "TDycore Mesh Examples\n\n";
2adcf2cb4SMatthew G. Knepley 
3adcf2cb4SMatthew G. Knepley #include <petscdmplex.h>
4adcf2cb4SMatthew G. Knepley 
5adcf2cb4SMatthew G. Knepley typedef struct {
6adcf2cb4SMatthew G. Knepley   PetscBool adapt; /* Flag for adaptation of the surface mesh */
7adcf2cb4SMatthew G. Knepley } AppCtx;
8adcf2cb4SMatthew G. Knepley 
9*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
10adcf2cb4SMatthew G. Knepley   PetscFunctionBeginUser;
11adcf2cb4SMatthew G. Knepley   options->adapt = PETSC_FALSE;
12adcf2cb4SMatthew G. Knepley 
13d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");
149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-adapt", "Flag for adaptation of the surface mesh", "ex10.c", options->adapt, &options->adapt, NULL));
15d0609cedSBarry Smith   PetscOptionsEnd();
16adcf2cb4SMatthew G. Knepley   PetscFunctionReturn(0);
17adcf2cb4SMatthew G. Knepley }
18adcf2cb4SMatthew G. Knepley 
19*9371c9d4SSatish Balay static PetscErrorCode CreateDomainLabel(DM dm) {
20adcf2cb4SMatthew G. Knepley   DMLabel  label;
21adcf2cb4SMatthew G. Knepley   PetscInt cStart, cEnd, c;
22adcf2cb4SMatthew G. Knepley 
23adcf2cb4SMatthew G. Knepley   PetscFunctionBeginUser;
248fb5bd83SMatthew G. Knepley   PetscCall(DMGetCoordinatesLocalSetUp(dm));
259566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "Cell Sets"));
269566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "Cell Sets", &label));
279566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
28adcf2cb4SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
29adcf2cb4SMatthew G. Knepley     PetscReal centroid[3], volume, x, y;
30adcf2cb4SMatthew G. Knepley 
319566063dSJacob Faibussowitsch     PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &volume, centroid, NULL));
32*9371c9d4SSatish Balay     x = centroid[0];
33*9371c9d4SSatish Balay     y = centroid[1];
34adcf2cb4SMatthew G. Knepley     /* Headwaters are (0.0,0.25)--(0.1,0.75) */
35*9371c9d4SSatish Balay     if ((x >= 0.0 && x < 0.1) && (y >= 0.25 && y <= 0.75)) {
36*9371c9d4SSatish Balay       PetscCall(DMLabelSetValue(label, c, 1));
37*9371c9d4SSatish Balay       continue;
38*9371c9d4SSatish Balay     }
39adcf2cb4SMatthew G. Knepley     /* River channel is (0.1,0.45)--(1.0,0.55) */
40*9371c9d4SSatish Balay     if ((x >= 0.1 && x <= 1.0) && (y >= 0.45 && y <= 0.55)) {
41*9371c9d4SSatish Balay       PetscCall(DMLabelSetValue(label, c, 2));
42*9371c9d4SSatish Balay       continue;
43*9371c9d4SSatish Balay     }
44adcf2cb4SMatthew G. Knepley   }
45adcf2cb4SMatthew G. Knepley   PetscFunctionReturn(0);
46adcf2cb4SMatthew G. Knepley }
47adcf2cb4SMatthew G. Knepley 
48*9371c9d4SSatish Balay static PetscErrorCode AdaptMesh(DM *dm, AppCtx *ctx) {
49adcf2cb4SMatthew G. Knepley   DM              dmCur = *dm;
50adcf2cb4SMatthew G. Knepley   DMLabel         label;
51adcf2cb4SMatthew G. Knepley   IS              valueIS, vIS;
52adcf2cb4SMatthew G. Knepley   PetscBool       hasLabel;
53adcf2cb4SMatthew G. Knepley   const PetscInt *values;
54adcf2cb4SMatthew G. Knepley   PetscReal      *volConst; /* Volume constraints for each label value */
55adcf2cb4SMatthew G. Knepley   PetscReal       ratio;
56adcf2cb4SMatthew G. Knepley   PetscInt        dim, Nv, v, cStart, cEnd, c;
57adcf2cb4SMatthew G. Knepley   PetscBool       adapt = PETSC_TRUE;
58adcf2cb4SMatthew G. Knepley 
59adcf2cb4SMatthew G. Knepley   PetscFunctionBeginUser;
60adcf2cb4SMatthew G. Knepley   if (!ctx->adapt) PetscFunctionReturn(0);
619566063dSJacob Faibussowitsch   PetscCall(DMHasLabel(*dm, "Cell Sets", &hasLabel));
629566063dSJacob Faibussowitsch   if (!hasLabel) PetscCall(CreateDomainLabel(*dm));
639566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(*dm, &dim));
64adcf2cb4SMatthew G. Knepley   ratio = PetscPowRealInt(0.5, dim);
65adcf2cb4SMatthew G. Knepley   /* Get volume constraints */
669566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(*dm, "Cell Sets", &label));
679566063dSJacob Faibussowitsch   PetscCall(DMLabelGetValueIS(label, &vIS));
689566063dSJacob Faibussowitsch   PetscCall(ISDuplicate(vIS, &valueIS));
699566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&vIS));
70adcf2cb4SMatthew G. Knepley   /* Sorting ruins the label */
719566063dSJacob Faibussowitsch   PetscCall(ISSort(valueIS));
729566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(valueIS, &Nv));
739566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(valueIS, &values));
749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nv, &volConst));
75adcf2cb4SMatthew G. Knepley   for (v = 0; v < Nv; ++v) {
76adcf2cb4SMatthew G. Knepley     char opt[128];
77adcf2cb4SMatthew G. Knepley 
78adcf2cb4SMatthew G. Knepley     volConst[v] = PETSC_MAX_REAL;
799566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(opt, 128, "-volume_constraint_%d", (int)values[v]));
809566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetReal(NULL, NULL, opt, &volConst[v], NULL));
81adcf2cb4SMatthew G. Knepley   }
829566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(valueIS, &values));
839566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&valueIS));
84adcf2cb4SMatthew G. Knepley   /* Adapt mesh iteratively */
85adcf2cb4SMatthew G. Knepley   while (adapt) {
86adcf2cb4SMatthew G. Knepley     DM       dmAdapt;
87adcf2cb4SMatthew G. Knepley     DMLabel  adaptLabel;
88adcf2cb4SMatthew G. Knepley     PetscInt nAdaptLoc[2], nAdapt[2];
89adcf2cb4SMatthew G. Knepley 
90adcf2cb4SMatthew G. Knepley     adapt        = PETSC_FALSE;
91adcf2cb4SMatthew G. Knepley     nAdaptLoc[0] = nAdaptLoc[1] = 0;
92adcf2cb4SMatthew G. Knepley     nAdapt[0] = nAdapt[1] = 0;
93adcf2cb4SMatthew G. Knepley     /* Adaptation is not preserving the domain label */
949566063dSJacob Faibussowitsch     PetscCall(DMHasLabel(dmCur, "Cell Sets", &hasLabel));
959566063dSJacob Faibussowitsch     if (!hasLabel) PetscCall(CreateDomainLabel(dmCur));
969566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dmCur, "Cell Sets", &label));
979566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValueIS(label, &vIS));
989566063dSJacob Faibussowitsch     PetscCall(ISDuplicate(vIS, &valueIS));
999566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&vIS));
100adcf2cb4SMatthew G. Knepley     /* Sorting directly the label's value IS would corrupt the label so we duplicate the IS first */
1019566063dSJacob Faibussowitsch     PetscCall(ISSort(valueIS));
1029566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(valueIS, &Nv));
1039566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(valueIS, &values));
104adcf2cb4SMatthew G. Knepley     /* Construct adaptation label */
1059566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
1069566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dmCur, 0, &cStart, &cEnd));
107adcf2cb4SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
108adcf2cb4SMatthew G. Knepley       PetscReal volume, centroid[3];
109adcf2cb4SMatthew G. Knepley       PetscInt  value, vidx;
110adcf2cb4SMatthew G. Knepley 
1119566063dSJacob Faibussowitsch       PetscCall(DMPlexComputeCellGeometryFVM(dmCur, c, &volume, centroid, NULL));
1129566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(label, c, &value));
113adcf2cb4SMatthew G. Knepley       if (value < 0) continue;
1149566063dSJacob Faibussowitsch       PetscCall(PetscFindInt(value, Nv, values, &vidx));
11563a3b9bcSJacob Faibussowitsch       PetscCheck(vidx >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %" PetscInt_FMT " for cell %" PetscInt_FMT " does not exist in label", value, c);
116*9371c9d4SSatish Balay       if (volume > volConst[vidx]) {
117*9371c9d4SSatish Balay         PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE));
118*9371c9d4SSatish Balay         ++nAdaptLoc[0];
119*9371c9d4SSatish Balay       }
120*9371c9d4SSatish Balay       if (volume < volConst[vidx] * ratio) {
121*9371c9d4SSatish Balay         PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_COARSEN));
122*9371c9d4SSatish Balay         ++nAdaptLoc[1];
123*9371c9d4SSatish Balay       }
124adcf2cb4SMatthew G. Knepley     }
1259566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(valueIS, &values));
1269566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&valueIS));
1279566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(&nAdaptLoc, &nAdapt, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dmCur)));
128adcf2cb4SMatthew G. Knepley     if (nAdapt[0]) {
12963a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(dmCur, "Adapted mesh, marking %" PetscInt_FMT " cells for refinement, and %" PetscInt_FMT " cells for coarsening\n", nAdapt[0], nAdapt[1]));
1309566063dSJacob Faibussowitsch       PetscCall(DMAdaptLabel(dmCur, adaptLabel, &dmAdapt));
1319566063dSJacob Faibussowitsch       PetscCall(DMDestroy(&dmCur));
1329566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(dmAdapt, NULL, "-adapt_dm_view"));
133adcf2cb4SMatthew G. Knepley       dmCur = dmAdapt;
134adcf2cb4SMatthew G. Knepley       adapt = PETSC_TRUE;
135adcf2cb4SMatthew G. Knepley     }
1369566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&adaptLabel));
137adcf2cb4SMatthew G. Knepley   }
1389566063dSJacob Faibussowitsch   PetscCall(PetscFree(volConst));
139adcf2cb4SMatthew G. Knepley   *dm = dmCur;
140adcf2cb4SMatthew G. Knepley   PetscFunctionReturn(0);
141adcf2cb4SMatthew G. Knepley }
142adcf2cb4SMatthew G. Knepley 
143*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
144adcf2cb4SMatthew G. Knepley   PetscInt dim;
145adcf2cb4SMatthew G. Knepley 
146adcf2cb4SMatthew G. Knepley   PetscFunctionBeginUser;
147adcf2cb4SMatthew G. Knepley   /* Create top surface */
1489566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
1499566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1509566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "init_"));
1519566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
1529566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, NULL));
153adcf2cb4SMatthew G. Knepley   /* Adapt surface */
1549566063dSJacob Faibussowitsch   PetscCall(AdaptMesh(dm, user));
155adcf2cb4SMatthew G. Knepley   /* Extrude surface to get volume mesh */
1569566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(*dm, &dim));
1579566063dSJacob Faibussowitsch   PetscCall(DMLocalizeCoordinates(*dm));
1589566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*dm, "Mesh"));
1599566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
1609566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
161adcf2cb4SMatthew G. Knepley   PetscFunctionReturn(0);
162adcf2cb4SMatthew G. Knepley }
163adcf2cb4SMatthew G. Knepley 
164*9371c9d4SSatish Balay int main(int argc, char **argv) {
165adcf2cb4SMatthew G. Knepley   DM     dm;
166adcf2cb4SMatthew G. Knepley   AppCtx user;
167adcf2cb4SMatthew G. Knepley 
168327415f7SBarry Smith   PetscFunctionBeginUser;
1699566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1709566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
1719566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1729566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
1739566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
174b122ec5aSJacob Faibussowitsch   return 0;
175adcf2cb4SMatthew G. Knepley }
176adcf2cb4SMatthew G. Knepley 
177adcf2cb4SMatthew G. Knepley /*TEST
178adcf2cb4SMatthew G. Knepley 
179adcf2cb4SMatthew G. Knepley   test:
180adcf2cb4SMatthew G. Knepley     suffix: 0
181adcf2cb4SMatthew G. Knepley     requires: triangle
182d410b0cfSMatthew G. Knepley     args: -init_dm_plex_dim 2 -init_dm_plex_box_faces 1,1 -dm_extrude 1 -dm_view
183adcf2cb4SMatthew G. Knepley 
184adcf2cb4SMatthew G. Knepley   test: # Regularly refine the surface before extrusion
185adcf2cb4SMatthew G. Knepley     suffix: 1
186adcf2cb4SMatthew G. Knepley     requires: triangle
187d410b0cfSMatthew G. Knepley     args: -init_dm_plex_dim 2 -init_dm_refine 2 -dm_extrude 1 -dm_view
188adcf2cb4SMatthew G. Knepley 
189adcf2cb4SMatthew G. Knepley   test: # Parallel run
190adcf2cb4SMatthew G. Knepley     suffix: 2
191adcf2cb4SMatthew G. Knepley     requires: triangle
192adcf2cb4SMatthew G. Knepley     nsize: 5
193e600fa54SMatthew G. Knepley     args: -init_dm_plex_dim 2 -init_dm_refine 3 -petscpartitioner_type simple -dm_extrude 3 -dm_view
194adcf2cb4SMatthew G. Knepley 
195adcf2cb4SMatthew G. Knepley   test: # adaptively refine the surface before extrusion
196adcf2cb4SMatthew G. Knepley     suffix: 3
197adcf2cb4SMatthew G. Knepley     requires: triangle
198d410b0cfSMatthew G. Knepley     args: -init_dm_plex_dim 2 -init_dm_plex_box_faces 5,5 -adapt -volume_constraint_1 0.01 -volume_constraint_2 0.000625 -dm_extrude 10
199adcf2cb4SMatthew G. Knepley 
200adcf2cb4SMatthew G. Knepley TEST*/
201