xref: /petsc/src/dm/impls/plex/tutorials/ex10.c (revision e600fa544e2bb197ca2af9b6e65ea465976dec56)
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 
9adcf2cb4SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
10adcf2cb4SMatthew G. Knepley {
11adcf2cb4SMatthew G. Knepley   PetscErrorCode ierr;
12adcf2cb4SMatthew G. Knepley 
13adcf2cb4SMatthew G. Knepley   PetscFunctionBeginUser;
14adcf2cb4SMatthew G. Knepley   options->adapt = PETSC_FALSE;
15adcf2cb4SMatthew G. Knepley 
16adcf2cb4SMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr);
17adcf2cb4SMatthew G. Knepley   ierr = PetscOptionsBool("-adapt", "Flag for adaptation of the surface mesh", "ex10.c", options->adapt, &options->adapt, NULL);CHKERRQ(ierr);
181e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
19adcf2cb4SMatthew G. Knepley   PetscFunctionReturn(0);
20adcf2cb4SMatthew G. Knepley }
21adcf2cb4SMatthew G. Knepley 
22adcf2cb4SMatthew G. Knepley static PetscErrorCode CreateDomainLabel(DM dm)
23adcf2cb4SMatthew G. Knepley {
24adcf2cb4SMatthew G. Knepley   DMLabel        label;
25adcf2cb4SMatthew G. Knepley   PetscInt       cStart, cEnd, c;
26adcf2cb4SMatthew G. Knepley   PetscErrorCode ierr;
27adcf2cb4SMatthew G. Knepley 
28adcf2cb4SMatthew G. Knepley   PetscFunctionBeginUser;
29adcf2cb4SMatthew G. Knepley   ierr = DMCreateLabel(dm, "Cell Sets");CHKERRQ(ierr);
30adcf2cb4SMatthew G. Knepley   ierr = DMGetLabel(dm, "Cell Sets", &label);CHKERRQ(ierr);
31adcf2cb4SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
32adcf2cb4SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
33adcf2cb4SMatthew G. Knepley     PetscReal centroid[3], volume, x, y;
34adcf2cb4SMatthew G. Knepley 
35adcf2cb4SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFVM(dm, c, &volume, centroid, NULL);CHKERRQ(ierr);
36adcf2cb4SMatthew G. Knepley     x = centroid[0]; y = centroid[1];
37adcf2cb4SMatthew G. Knepley     /* Headwaters are (0.0,0.25)--(0.1,0.75) */
38adcf2cb4SMatthew G. Knepley     if ((x >= 0.0 && x <  0.1) && (y >= 0.25 && y <= 0.75)) {ierr = DMLabelSetValue(label, c, 1);CHKERRQ(ierr);continue;}
39adcf2cb4SMatthew G. Knepley     /* River channel is (0.1,0.45)--(1.0,0.55) */
40adcf2cb4SMatthew G. Knepley     if ((x >= 0.1 && x <= 1.0) && (y >= 0.45 && y <= 0.55)) {ierr = DMLabelSetValue(label, c, 2);CHKERRQ(ierr);continue;}
41adcf2cb4SMatthew G. Knepley   }
42adcf2cb4SMatthew G. Knepley   PetscFunctionReturn(0);
43adcf2cb4SMatthew G. Knepley }
44adcf2cb4SMatthew G. Knepley 
45adcf2cb4SMatthew G. Knepley static PetscErrorCode AdaptMesh(DM *dm, AppCtx *ctx)
46adcf2cb4SMatthew G. Knepley {
47adcf2cb4SMatthew G. Knepley   DM              dmCur = *dm;
48adcf2cb4SMatthew G. Knepley   DMLabel         label;
49adcf2cb4SMatthew G. Knepley   IS              valueIS, vIS;
50adcf2cb4SMatthew G. Knepley   PetscBool       hasLabel;
51adcf2cb4SMatthew G. Knepley   const PetscInt *values;
52adcf2cb4SMatthew G. Knepley   PetscReal      *volConst; /* Volume constraints for each label value */
53adcf2cb4SMatthew G. Knepley   PetscReal       ratio;
54adcf2cb4SMatthew G. Knepley   PetscInt        dim, Nv, v, cStart, cEnd, c;
55adcf2cb4SMatthew G. Knepley   PetscBool       adapt = PETSC_TRUE;
56adcf2cb4SMatthew G. Knepley   PetscErrorCode  ierr;
57adcf2cb4SMatthew G. Knepley 
58adcf2cb4SMatthew G. Knepley   PetscFunctionBeginUser;
59adcf2cb4SMatthew G. Knepley   if (!ctx->adapt) PetscFunctionReturn(0);
60adcf2cb4SMatthew G. Knepley   ierr = DMHasLabel(*dm, "Cell Sets", &hasLabel);CHKERRQ(ierr);
61adcf2cb4SMatthew G. Knepley   if (!hasLabel) {ierr = CreateDomainLabel(*dm);CHKERRQ(ierr);}
62adcf2cb4SMatthew G. Knepley   ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr);
63adcf2cb4SMatthew G. Knepley   ratio = PetscPowRealInt(0.5, dim);
64adcf2cb4SMatthew G. Knepley   /* Get volume constraints */
65adcf2cb4SMatthew G. Knepley   ierr = DMGetLabel(*dm, "Cell Sets", &label);CHKERRQ(ierr);
66adcf2cb4SMatthew G. Knepley   ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
67adcf2cb4SMatthew G. Knepley   ierr = ISDuplicate(vIS, &valueIS);CHKERRQ(ierr);
68adcf2cb4SMatthew G. Knepley   ierr = ISDestroy(&vIS);CHKERRQ(ierr);
69adcf2cb4SMatthew G. Knepley   /* Sorting ruins the label */
70adcf2cb4SMatthew G. Knepley   ierr = ISSort(valueIS);CHKERRQ(ierr);
71adcf2cb4SMatthew G. Knepley   ierr = ISGetLocalSize(valueIS, &Nv);CHKERRQ(ierr);
72adcf2cb4SMatthew G. Knepley   ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr);
73adcf2cb4SMatthew G. Knepley   ierr = PetscMalloc1(Nv, &volConst);CHKERRQ(ierr);
74adcf2cb4SMatthew G. Knepley   for (v = 0; v < Nv; ++v) {
75adcf2cb4SMatthew G. Knepley     char opt[128];
76adcf2cb4SMatthew G. Knepley 
77adcf2cb4SMatthew G. Knepley     volConst[v] = PETSC_MAX_REAL;
78adcf2cb4SMatthew G. Knepley     ierr = PetscSNPrintf(opt, 128, "-volume_constraint_%d", (int) values[v]);CHKERRQ(ierr);
79adcf2cb4SMatthew G. Knepley     ierr = PetscOptionsGetReal(NULL, NULL, opt, &volConst[v], NULL);CHKERRQ(ierr);
80adcf2cb4SMatthew G. Knepley   }
81adcf2cb4SMatthew G. Knepley   ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr);
82adcf2cb4SMatthew G. Knepley   ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
83adcf2cb4SMatthew G. Knepley   /* Adapt mesh iteratively */
84adcf2cb4SMatthew G. Knepley   while (adapt) {
85adcf2cb4SMatthew G. Knepley     DM       dmAdapt;
86adcf2cb4SMatthew G. Knepley     DMLabel  adaptLabel;
87adcf2cb4SMatthew G. Knepley     PetscInt nAdaptLoc[2], nAdapt[2];
88adcf2cb4SMatthew G. Knepley 
89adcf2cb4SMatthew G. Knepley     adapt = PETSC_FALSE;
90adcf2cb4SMatthew G. Knepley     nAdaptLoc[0] = nAdaptLoc[1] = 0;
91adcf2cb4SMatthew G. Knepley     nAdapt[0]    = nAdapt[1]    = 0;
92adcf2cb4SMatthew G. Knepley     /* Adaptation is not preserving the domain label */
93adcf2cb4SMatthew G. Knepley     ierr = DMHasLabel(dmCur, "Cell Sets", &hasLabel);CHKERRQ(ierr);
94adcf2cb4SMatthew G. Knepley     if (!hasLabel) {ierr = CreateDomainLabel(dmCur);CHKERRQ(ierr);}
95adcf2cb4SMatthew G. Knepley     ierr = DMGetLabel(dmCur, "Cell Sets", &label);CHKERRQ(ierr);
96adcf2cb4SMatthew G. Knepley     ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr);
97adcf2cb4SMatthew G. Knepley     ierr = ISDuplicate(vIS, &valueIS);CHKERRQ(ierr);
98adcf2cb4SMatthew G. Knepley     ierr = ISDestroy(&vIS);CHKERRQ(ierr);
99adcf2cb4SMatthew G. Knepley     /* Sorting directly the label's value IS would corrupt the label so we duplicate the IS first */
100adcf2cb4SMatthew G. Knepley     ierr = ISSort(valueIS);CHKERRQ(ierr);
101adcf2cb4SMatthew G. Knepley     ierr = ISGetLocalSize(valueIS, &Nv);CHKERRQ(ierr);
102adcf2cb4SMatthew G. Knepley     ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr);
103adcf2cb4SMatthew G. Knepley     /* Construct adaptation label */
104adcf2cb4SMatthew G. Knepley     ierr = DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel);CHKERRQ(ierr);
105adcf2cb4SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dmCur, 0, &cStart, &cEnd);CHKERRQ(ierr);
106adcf2cb4SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
107adcf2cb4SMatthew G. Knepley       PetscReal volume, centroid[3];
108adcf2cb4SMatthew G. Knepley       PetscInt  value, vidx;
109adcf2cb4SMatthew G. Knepley 
110adcf2cb4SMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dmCur, c, &volume, centroid, NULL);CHKERRQ(ierr);
111adcf2cb4SMatthew G. Knepley       ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr);
112adcf2cb4SMatthew G. Knepley       if (value < 0) continue;
113adcf2cb4SMatthew G. Knepley       ierr = PetscFindInt(value, Nv, values, &vidx);CHKERRQ(ierr);
1142c71b3e2SJacob Faibussowitsch       PetscCheckFalse(vidx < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %D for cell %D does not exist in label", value, c);
115adcf2cb4SMatthew G. Knepley       if (volume > volConst[vidx])        {ierr = DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE);CHKERRQ(ierr);  ++nAdaptLoc[0];}
116adcf2cb4SMatthew G. Knepley       if (volume < volConst[vidx]*ratio) {ierr = DMLabelSetValue(adaptLabel, c, DM_ADAPT_COARSEN);CHKERRQ(ierr); ++nAdaptLoc[1];}
117adcf2cb4SMatthew G. Knepley     }
118adcf2cb4SMatthew G. Knepley     ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr);
119adcf2cb4SMatthew G. Knepley     ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
12055b25c41SPierre Jolivet     ierr = MPI_Allreduce(&nAdaptLoc, &nAdapt, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) dmCur));CHKERRMPI(ierr);
121adcf2cb4SMatthew G. Knepley     if (nAdapt[0]) {
1227d3de750SJacob Faibussowitsch       ierr = PetscInfo(dmCur, "Adapted mesh, marking %D cells for refinement, and %D cells for coarsening\n", nAdapt[0], nAdapt[1]);CHKERRQ(ierr);
123adcf2cb4SMatthew G. Knepley       ierr = DMAdaptLabel(dmCur, adaptLabel, &dmAdapt);CHKERRQ(ierr);
1241e1ea65dSPierre Jolivet       ierr = DMDestroy(&dmCur);CHKERRQ(ierr);
125adcf2cb4SMatthew G. Knepley       ierr = DMViewFromOptions(dmAdapt, NULL, "-adapt_dm_view");CHKERRQ(ierr);
126adcf2cb4SMatthew G. Knepley       dmCur = dmAdapt;
127adcf2cb4SMatthew G. Knepley       adapt = PETSC_TRUE;
128adcf2cb4SMatthew G. Knepley     }
129adcf2cb4SMatthew G. Knepley     ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr);
130adcf2cb4SMatthew G. Knepley   }
131adcf2cb4SMatthew G. Knepley   ierr = PetscFree(volConst);CHKERRQ(ierr);
132adcf2cb4SMatthew G. Knepley   *dm = dmCur;
133adcf2cb4SMatthew G. Knepley   PetscFunctionReturn(0);
134adcf2cb4SMatthew G. Knepley }
135adcf2cb4SMatthew G. Knepley 
136adcf2cb4SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
137adcf2cb4SMatthew G. Knepley {
138adcf2cb4SMatthew G. Knepley   PetscInt       dim;
139adcf2cb4SMatthew G. Knepley   PetscErrorCode ierr;
140adcf2cb4SMatthew G. Knepley 
141adcf2cb4SMatthew G. Knepley   PetscFunctionBeginUser;
142adcf2cb4SMatthew G. Knepley   /* Create top surface */
14330602db0SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
14430602db0SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
14530602db0SMatthew G. Knepley   ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, "init_");CHKERRQ(ierr);
14630602db0SMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
147d410b0cfSMatthew G. Knepley   ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, NULL);CHKERRQ(ierr);
148adcf2cb4SMatthew G. Knepley   /* Adapt surface */
149adcf2cb4SMatthew G. Knepley   ierr = AdaptMesh(dm, user);CHKERRQ(ierr);
150adcf2cb4SMatthew G. Knepley   /* Extrude surface to get volume mesh */
151adcf2cb4SMatthew G. Knepley   ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr);
152adcf2cb4SMatthew G. Knepley   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
153adcf2cb4SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
154adcf2cb4SMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
155adcf2cb4SMatthew G. Knepley   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
156adcf2cb4SMatthew G. Knepley   PetscFunctionReturn(0);
157adcf2cb4SMatthew G. Knepley }
158adcf2cb4SMatthew G. Knepley 
159adcf2cb4SMatthew G. Knepley int main(int argc, char **argv)
160adcf2cb4SMatthew G. Knepley {
161adcf2cb4SMatthew G. Knepley   DM             dm;
162adcf2cb4SMatthew G. Knepley   AppCtx         user;
163adcf2cb4SMatthew G. Knepley   PetscErrorCode ierr;
164adcf2cb4SMatthew G. Knepley 
165adcf2cb4SMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
166adcf2cb4SMatthew G. Knepley   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
167adcf2cb4SMatthew G. Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
168adcf2cb4SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
169adcf2cb4SMatthew G. Knepley   ierr = PetscFinalize();
170adcf2cb4SMatthew G. Knepley   return ierr;
171adcf2cb4SMatthew G. Knepley }
172adcf2cb4SMatthew G. Knepley 
173adcf2cb4SMatthew G. Knepley /*TEST
174adcf2cb4SMatthew G. Knepley 
175adcf2cb4SMatthew G. Knepley   test:
176adcf2cb4SMatthew G. Knepley     suffix: 0
177adcf2cb4SMatthew G. Knepley     requires: triangle
178d410b0cfSMatthew G. Knepley     args: -init_dm_plex_dim 2 -init_dm_plex_box_faces 1,1 -dm_extrude 1 -dm_view
179adcf2cb4SMatthew G. Knepley 
180adcf2cb4SMatthew G. Knepley   test: # Regularly refine the surface before extrusion
181adcf2cb4SMatthew G. Knepley     suffix: 1
182adcf2cb4SMatthew G. Knepley     requires: triangle
183d410b0cfSMatthew G. Knepley     args: -init_dm_plex_dim 2 -init_dm_refine 2 -dm_extrude 1 -dm_view
184adcf2cb4SMatthew G. Knepley 
185adcf2cb4SMatthew G. Knepley   test: # Parallel run
186adcf2cb4SMatthew G. Knepley     suffix: 2
187adcf2cb4SMatthew G. Knepley     requires: triangle
188adcf2cb4SMatthew G. Knepley     nsize: 5
189*e600fa54SMatthew G. Knepley     args: -init_dm_plex_dim 2 -init_dm_refine 3 -petscpartitioner_type simple -dm_extrude 3 -dm_view
190adcf2cb4SMatthew G. Knepley 
191adcf2cb4SMatthew G. Knepley   test: # adaptively refine the surface before extrusion
192adcf2cb4SMatthew G. Knepley     suffix: 3
193adcf2cb4SMatthew G. Knepley     requires: triangle
194d410b0cfSMatthew 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
195adcf2cb4SMatthew G. Knepley 
196adcf2cb4SMatthew G. Knepley TEST*/
197