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