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