xref: /petsc/src/dm/impls/plex/tests/ex60.c (revision c3e4dd79e17d3f0a51a524340fae4661630fd09e)
1 static char help[] = "Test metric utils in the uniform, isotropic case.\n\n";
2 
3 #include <petscdmplex.h>
4 
5 static PetscErrorCode bowl(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
6 {
7   PetscInt d;
8 
9   *u = 0.0;
10   for (d = 0; d < dim; d++) *u += 0.5*(x[d] - 0.5)*(x[d] - 0.5);
11 
12   return 0;
13 }
14 
15 static PetscErrorCode CreateIndicator(DM dm, Vec *indicator, DM *dmIndi)
16 {
17   MPI_Comm       comm;
18   PetscFE        fe;
19   PetscInt       dim;
20 
21   PetscFunctionBeginUser;
22   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
23   PetscCall(DMClone(dm, dmIndi));
24   PetscCall(DMGetDimension(dm, &dim));
25   PetscCall(PetscFECreateLagrange(comm, dim, 1, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
26   PetscCall(DMSetField(*dmIndi, 0, NULL, (PetscObject)fe));
27   PetscCall(DMCreateDS(*dmIndi));
28   PetscCall(PetscFEDestroy(&fe));
29   PetscCall(DMCreateLocalVector(*dmIndi, indicator));
30   PetscFunctionReturn(0);
31 }
32 
33 int main(int argc, char **argv) {
34   DM              dm, dmAdapt;
35   DMLabel         bdLabel = NULL, rgLabel = NULL;
36   MPI_Comm        comm;
37   PetscBool       uniform = PETSC_FALSE, isotropic = PETSC_FALSE, noTagging = PETSC_FALSE;
38   PetscInt        dim;
39   PetscReal       scaling = 1.0;
40   Vec             metric;
41 
42   /* Set up */
43   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
44   comm = PETSC_COMM_WORLD;
45   PetscOptionsBegin(comm, "", "Mesh adaptation options", "DMPLEX");
46   PetscCall(PetscOptionsBool("-noTagging", "Should tag preservation testing be turned off?", "ex60.c", noTagging, &noTagging, NULL));
47   PetscOptionsEnd();
48 
49   /* Create box mesh */
50   PetscCall(DMCreate(comm, &dm));
51   PetscCall(DMSetType(dm, DMPLEX));
52   PetscCall(DMSetFromOptions(dm));
53   PetscCall(PetscObjectSetName((PetscObject) dm, "DM_init"));
54   PetscCall(DMViewFromOptions(dm, NULL, "-initial_mesh_view"));
55   PetscCall(DMGetDimension(dm, &dim));
56 
57   /* Set tags to be preserved */
58   if (!noTagging) {
59     DM                 cdm;
60     PetscInt           cStart, cEnd, c, fStart, fEnd, f, vStart, vEnd;
61     const PetscScalar *coords;
62     Vec                coordinates;
63 
64     /* Cell tags */
65     PetscCall(DMCreateLabel(dm, "Cell Sets"));
66     PetscCall(DMGetLabel(dm, "Cell Sets", &rgLabel));
67     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
68     for (c = cStart; c < cEnd; ++c) {
69       PetscReal centroid[3], volume, x;
70 
71       PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &volume, centroid, NULL));
72       x = centroid[0];
73       if (x < 0.5) PetscCall(DMLabelSetValue(rgLabel, c, 3));
74       else         PetscCall(DMLabelSetValue(rgLabel, c, 4));
75     }
76 
77     /* Face tags */
78     PetscCall(DMCreateLabel(dm, "Face Sets"));
79     PetscCall(DMGetLabel(dm, "Face Sets", &bdLabel));
80     PetscCall(DMPlexMarkBoundaryFaces(dm, 1, bdLabel));
81     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
82     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
83     PetscCall(DMGetCoordinateDM(dm, &cdm));
84     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
85     PetscCall(VecGetArrayRead(coordinates, &coords));
86     for (f = fStart; f < fEnd; ++f) {
87       PetscBool flg = PETSC_TRUE;
88       PetscInt *closure = NULL, closureSize, cl;
89       PetscReal eps = 1.0e-08;
90 
91       PetscCall(DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
92       for (cl = 0; cl < closureSize*2; cl += 2) {
93         PetscInt   off = closure[cl];
94         PetscReal *x;
95 
96         if ((off < vStart) || (off >= vEnd)) continue;
97         PetscCall(DMPlexPointLocalRead(cdm, off, coords, &x));
98         if ((x[0] < 0.5 - eps) || (x[0] > 0.5 + eps)) flg = PETSC_FALSE;
99       }
100       if (flg) PetscCall(DMLabelSetValue(bdLabel, f, 2));
101       PetscCall(DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure));
102     }
103     PetscCall(VecRestoreArrayRead(coordinates, &coords));
104   }
105 
106   /* Construct metric */
107   PetscCall(DMPlexMetricSetFromOptions(dm));
108   PetscCall(DMPlexMetricIsUniform(dm, &uniform));
109   PetscCall(DMPlexMetricIsIsotropic(dm, &isotropic));
110   if (uniform) {
111     PetscCall(DMPlexMetricCreateUniform(dm, 0, scaling, &metric));
112   }
113   else {
114     DM  dmIndi;
115     Vec indicator;
116 
117     /* Construct "error indicator" */
118     PetscCall(CreateIndicator(dm, &indicator, &dmIndi));
119     if (isotropic) {
120 
121       /* Isotropic case: just specify unity */
122       PetscCall(VecSet(indicator, scaling));
123       PetscCall(DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric));
124 
125     } else {
126       PetscFE fe;
127 
128       /* 'Anisotropic' case: approximate the identity by recovering the Hessian of a parabola */
129       DM               dmGrad;
130       PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar*, void*) = {bowl};
131       Vec              gradient;
132 
133       /* Project the parabola into P1 space */
134       PetscCall(DMProjectFunctionLocal(dmIndi, 0.0, funcs, NULL, INSERT_ALL_VALUES, indicator));
135 
136       /* Approximate the gradient */
137       PetscCall(DMClone(dmIndi, &dmGrad));
138       PetscCall(PetscFECreateLagrange(comm, dim, dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe));
139       PetscCall(DMSetField(dmGrad, 0, NULL, (PetscObject)fe));
140       PetscCall(DMCreateDS(dmGrad));
141       PetscCall(PetscFEDestroy(&fe));
142       PetscCall(DMCreateLocalVector(dmGrad, &gradient));
143       PetscCall(DMPlexComputeGradientClementInterpolant(dmIndi, indicator, gradient));
144       PetscCall(VecViewFromOptions(gradient, NULL, "-adapt_gradient_view"));
145 
146       /* Approximate the Hessian */
147       PetscCall(DMPlexMetricCreate(dm, 0, &metric));
148       PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, gradient, metric));
149       PetscCall(VecViewFromOptions(metric, NULL, "-adapt_hessian_view"));
150       PetscCall(VecDestroy(&gradient));
151       PetscCall(DMDestroy(&dmGrad));
152     }
153     PetscCall(VecDestroy(&indicator));
154     PetscCall(DMDestroy(&dmIndi));
155   }
156 
157   /* Test metric routines */
158   {
159     DM        dmDet;
160     PetscReal errornorm, norm, tol = 1.0e-10, weights[2] = {0.8, 0.2};
161     Vec       metric1, metric2, metricComb, determinant;
162     Vec       metrics[2];
163 
164     PetscCall(VecDuplicate(metric, &metric1));
165     PetscCall(VecSet(metric1, 0));
166     PetscCall(VecAXPY(metric1, 0.625, metric));
167     PetscCall(VecDuplicate(metric, &metric2));
168     PetscCall(VecSet(metric2, 0));
169     PetscCall(VecAXPY(metric2, 2.5, metric));
170     metrics[0] = metric1;
171     metrics[1] = metric2;
172 
173     /* Test metric average */
174     PetscCall(DMPlexMetricCreate(dm, 0, &metricComb));
175     PetscCall(DMPlexMetricAverage(dm, 2, weights, metrics, metricComb));
176     PetscCall(VecAXPY(metricComb, -1, metric));
177     PetscCall(VecNorm(metric, NORM_2, &norm));
178     PetscCall(VecNorm(metricComb, NORM_2, &errornorm));
179     errornorm /= norm;
180     PetscCall(PetscPrintf(comm, "Metric average L2 error: %.4f%%\n", (double)(100*errornorm)));
181     PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed");
182 
183     /* Test metric intersection */
184     if (isotropic) {
185       PetscCall(DMPlexMetricIntersection(dm, 2, metrics, metricComb));
186       PetscCall(VecAXPY(metricComb, -1, metric2));
187       PetscCall(VecNorm(metricComb, NORM_2, &errornorm));
188       errornorm /= norm;
189       PetscCall(PetscPrintf(comm, "Metric intersection L2 error: %.4f%%\n", (double)(100*errornorm)));
190       PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed");
191     }
192     PetscCall(VecDestroy(&metric2));
193     PetscCall(VecDestroy(&metricComb));
194 
195     /* Test metric SPD enforcement */
196     PetscCall(DMPlexMetricDeterminantCreate(dm, 0, &determinant, &dmDet));
197     PetscCall(DMPlexMetricEnforceSPD(dm, metric, PETSC_TRUE, PETSC_TRUE, metric1, determinant));
198     if (isotropic) {
199       Vec err;
200 
201       PetscCall(VecDuplicate(determinant, &err));
202       PetscCall(VecSet(err, 1.0));
203       PetscCall(VecNorm(err, NORM_2, &norm));
204       PetscCall(VecAXPY(err, -1, determinant));
205       PetscCall(VecNorm(err, NORM_2, &errornorm));
206       PetscCall(VecDestroy(&err));
207       errornorm /= norm;
208       PetscCall(PetscPrintf(comm, "Metric determinant L2 error: %.4f%%\n", (double)(100*errornorm)));
209       PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Determinant is not unit");
210       PetscCall(VecAXPY(metric1, -1, metric));
211       PetscCall(VecNorm(metric1, NORM_2, &errornorm));
212       errornorm /= norm;
213       PetscCall(PetscPrintf(comm, "Metric SPD enforcement L2 error: %.4f%%\n", (double)(100*errornorm)));
214       PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed");
215     }
216 
217     /* Test metric normalization */
218     PetscCall(DMPlexMetricNormalize(dm, metric, PETSC_TRUE, PETSC_TRUE, metric1, determinant));
219     if (isotropic) {
220       PetscReal target;
221 
222       PetscCall(DMPlexMetricGetTargetComplexity(dm, &target));
223       scaling = PetscPowReal(target, 2.0/dim);
224       if (uniform) {
225         PetscCall(DMPlexMetricCreateUniform(dm, 0, scaling, &metric2));
226       } else {
227         DM  dmIndi;
228         Vec indicator;
229 
230         PetscCall(CreateIndicator(dm, &indicator, &dmIndi));
231         PetscCall(VecSet(indicator, scaling));
232         PetscCall(DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric2));
233         PetscCall(DMDestroy(&dmIndi));
234         PetscCall(VecDestroy(&indicator));
235       }
236       PetscCall(VecAXPY(metric2, -1, metric1));
237       PetscCall(VecNorm(metric2, NORM_2, &errornorm));
238       errornorm /= norm;
239       PetscCall(PetscPrintf(comm, "Metric normalization L2 error: %.4f%%\n", (double)(100*errornorm)));
240       PetscCheck(errornorm < tol, comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed");
241     }
242     PetscCall(VecDestroy(&determinant));
243     PetscCall(DMDestroy(&dmDet));
244     PetscCall(VecCopy(metric1, metric));
245     PetscCall(VecDestroy(&metric2));
246     PetscCall(VecDestroy(&metric1));
247   }
248 
249   /* Adapt the mesh */
250   PetscCall(DMAdaptMetric(dm, metric, bdLabel, rgLabel, &dmAdapt));
251   PetscCall(DMDestroy(&dm));
252   PetscCall(PetscObjectSetName((PetscObject) dmAdapt, "DM_adapted"));
253   PetscCall(VecDestroy(&metric));
254   PetscCall(DMViewFromOptions(dmAdapt, NULL, "-adapted_mesh_view"));
255 
256   /* Test tag preservation */
257   if (!noTagging) {
258     PetscBool hasTag;
259     PetscInt  size;
260 
261     PetscCall(DMGetLabel(dmAdapt, "Face Sets", &bdLabel));
262     PetscCall(DMLabelHasStratum(bdLabel, 1, &hasTag));
263     PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 1");
264     PetscCall(DMLabelHasStratum(bdLabel, 2, &hasTag));
265     PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have face tag 2");
266     PetscCall(DMLabelGetNumValues(bdLabel, &size));
267     PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of face tags (got %" PetscInt_FMT ", expected 2)", size);
268 
269     PetscCall(DMGetLabel(dmAdapt, "Cell Sets", &rgLabel));
270     PetscCall(DMLabelHasStratum(rgLabel, 3, &hasTag));
271     PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 3");
272     PetscCall(DMLabelHasStratum(rgLabel, 4, &hasTag));
273     PetscCheck(hasTag, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh does not have cell tag 4");
274     PetscCall(DMLabelGetNumValues(rgLabel, &size));
275     PetscCheck(size == 2, comm, PETSC_ERR_ARG_OUTOFRANGE, "Adapted mesh has the wrong number of cell tags (got %" PetscInt_FMT ", expected 2)", size);
276   }
277 
278   /* Clean up */
279   PetscCall(DMDestroy(&dmAdapt));
280   PetscCall(PetscFinalize());
281   return 0;
282 }
283 
284 /*TEST
285 
286   testset:
287     requires: pragmatic
288     args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging
289 
290     test:
291       suffix: uniform_2d_pragmatic
292       args: -dm_plex_metric_uniform
293     test:
294       suffix: iso_2d_pragmatic
295       args: -dm_plex_metric_isotropic
296     test:
297       suffix: hessian_2d_pragmatic
298 
299   testset:
300     requires: pragmatic tetgen
301     args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor pragmatic -noTagging
302 
303     test:
304       suffix: uniform_3d_pragmatic
305       args: -dm_plex_metric_uniform -noTagging
306     test:
307       suffix: iso_3d_pragmatic
308       args: -dm_plex_metric_isotropic -noTagging
309     test:
310       suffix: hessian_3d_pragmatic
311 
312   testset:
313     requires: mmg
314     args: -dm_plex_box_faces 4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg
315 
316     test:
317       suffix: uniform_2d_mmg
318       args: -dm_plex_metric_uniform
319     test:
320       suffix: iso_2d_mmg
321       args: -dm_plex_metric_isotropic
322     test:
323       suffix: hessian_2d_mmg
324 
325   testset:
326     requires: mmg tetgen
327     args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor mmg
328 
329     test:
330       suffix: uniform_3d_mmg
331       args: -dm_plex_metric_uniform
332     test:
333       suffix: iso_3d_mmg
334       args: -dm_plex_metric_isotropic
335     test:
336       suffix: hessian_3d_mmg
337 
338   testset:
339     requires: parmmg tetgen
340     nsize: 2
341     args: -dm_plex_dim 3 -dm_plex_box_faces 4,4,4 -dm_plex_metric_target_complexity 100 -dm_adaptor parmmg
342 
343     test:
344       suffix: uniform_3d_parmmg
345       args: -dm_plex_metric_uniform
346     test:
347       suffix: iso_3d_parmmg
348       args: -dm_plex_metric_isotropic
349     test:
350       suffix: hessian_3d_parmmg
351 
352 TEST*/
353