xref: /petsc/src/dm/impls/plex/tests/ex60.c (revision e831869d00a6b6bf5bd64f1ce9662fca56500707)
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 int main(int argc, char **argv) {
16   DM              dm, dmDist, dmAdapt;
17   DMLabel         bdLabel = NULL;
18   MPI_Comm        comm;
19   PetscBool       uniform = PETSC_FALSE, isotropic = PETSC_FALSE;
20   PetscErrorCode  ierr;
21   PetscInt       *faces, dim = 3, numEdges = 4, d;
22   PetscReal       scaling = 1.0;
23   Vec             metric;
24 
25   /* Set up */
26   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
27   comm = PETSC_COMM_WORLD;
28   ierr = PetscOptionsBegin(comm, "", "Mesh adaptation options", "DMPLEX");CHKERRQ(ierr);
29   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex60.c", dim, &dim, NULL, 2, 3);CHKERRQ(ierr);
30   ierr = PetscOptionsBoundedInt("-num_edges", "Number of edges on each boundary of the initial mesh", "ex60.c", numEdges, &numEdges, NULL, 0);CHKERRQ(ierr);
31   ierr = PetscOptionsBool("-uniform", "Should the metric be assumed uniform?", "ex60.c", uniform, &uniform, NULL);CHKERRQ(ierr);
32   ierr = PetscOptionsBool("-isotropic", "Should the metric be assumed isotropic, or computed as a recovered Hessian?", "ex60.c", isotropic, &isotropic, NULL);CHKERRQ(ierr);
33   ierr = PetscOptionsEnd();
34 
35   /* Create box mesh */
36   ierr = PetscMalloc1(dim, &faces);CHKERRQ(ierr);
37   for (d = 0; d < dim; ++d) faces[d] = numEdges;
38   ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_TRUE, faces, NULL, NULL, NULL, PETSC_TRUE, &dm);CHKERRQ(ierr);
39   ierr = PetscFree(faces);CHKERRQ(ierr);
40   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
41 
42   /* Distribute mesh over processes */
43   ierr = DMPlexDistribute(dm, 0, NULL, &dmDist);CHKERRQ(ierr);
44   if (dmDist) {
45     ierr = DMDestroy(&dm);CHKERRQ(ierr);
46     dm = dmDist;
47   }
48   ierr = PetscObjectSetName((PetscObject) dm, "DM_init");CHKERRQ(ierr);
49   ierr = DMViewFromOptions(dm, NULL, "-initial_mesh_view");CHKERRQ(ierr);
50 
51   /* Construct metric */
52   if (uniform) {
53     if (isotropic) { ierr = DMPlexMetricCreateUniform(dm, 0, scaling, &metric);CHKERRQ(ierr); }
54     else SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Uniform anisotropic metrics not supported.");
55   }
56   else {
57     DM      dmIndi;
58     PetscFE fe;
59     Vec     indicator;
60 
61     /* Construct "error indicator" */
62     ierr = DMClone(dm, &dmIndi);CHKERRQ(ierr);
63     ierr = PetscFECreateLagrange(comm, dim, 1, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
64     ierr = DMSetField(dmIndi, 0, NULL, (PetscObject)fe);CHKERRQ(ierr);
65     ierr = DMCreateDS(dmIndi);CHKERRQ(ierr);
66     ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
67     ierr = DMCreateLocalVector(dmIndi, &indicator);CHKERRQ(ierr);
68     if (isotropic) {
69 
70       /* Isotropic case: just specify unity */
71       ierr = VecSet(indicator, scaling);CHKERRQ(ierr);
72       ierr = DMPlexMetricCreateIsotropic(dm, 0, indicator, &metric);CHKERRQ(ierr);
73 
74     } else {
75 
76       /* 'Anisotropic' case: approximate the identity by recovering the Hessian of a parabola */
77       DM               dmGrad;
78       PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar*, void*) = {bowl};
79       Vec              gradient;
80 
81       /* Project the parabola into P1 space */
82       ierr = DMProjectFunctionLocal(dmIndi, 0.0, funcs, NULL, INSERT_ALL_VALUES, indicator);CHKERRQ(ierr);
83 
84       /* Approximate the gradient */
85       ierr = DMClone(dmIndi, &dmGrad);CHKERRQ(ierr);
86       ierr = PetscFECreateLagrange(comm, dim, dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
87       ierr = DMSetField(dmGrad, 0, NULL, (PetscObject)fe);CHKERRQ(ierr);
88       ierr = DMCreateDS(dmGrad);CHKERRQ(ierr);
89       ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
90       ierr = DMCreateLocalVector(dmGrad, &gradient);CHKERRQ(ierr);
91       ierr = DMPlexComputeGradientClementInterpolant(dmIndi, indicator, gradient);CHKERRQ(ierr);
92       ierr = VecViewFromOptions(gradient, NULL, "-adapt_gradient_view");CHKERRQ(ierr);
93 
94       /* Approximate the Hessian */
95       ierr = DMPlexMetricCreate(dm, 0, &metric);CHKERRQ(ierr);
96       ierr = DMPlexComputeGradientClementInterpolant(dmGrad, gradient, metric);CHKERRQ(ierr);
97       ierr = VecViewFromOptions(metric, NULL, "-adapt_hessian_view");CHKERRQ(ierr);
98       ierr = VecDestroy(&gradient);CHKERRQ(ierr);
99       ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
100     }
101     ierr = VecDestroy(&indicator);CHKERRQ(ierr);
102     ierr = DMDestroy(&dmIndi);CHKERRQ(ierr);
103   }
104 
105   /* Test metric routines */
106   {
107     PetscReal errornorm, norm, tol = 1.0e-10, weights[2] = {0.8, 0.2};
108     Vec       metric1, metric2, metricComb;
109     Vec       metrics[2];
110 
111     ierr = VecDuplicate(metric, &metric1);CHKERRQ(ierr);
112     ierr = VecSet(metric1, 0);CHKERRQ(ierr);
113     ierr = VecAXPY(metric1, 0.625, metric);CHKERRQ(ierr);
114     ierr = VecDuplicate(metric, &metric2);CHKERRQ(ierr);
115     ierr = VecSet(metric2, 0);CHKERRQ(ierr);
116     ierr = VecAXPY(metric2, 2.5, metric);CHKERRQ(ierr);
117     metrics[0] = metric1;
118     metrics[1] = metric2;
119 
120     /* Test metric average */
121     ierr = DMPlexMetricAverage(dm, 2, weights, metrics, &metricComb);CHKERRQ(ierr);
122     ierr = VecAXPY(metricComb, -1, metric);CHKERRQ(ierr);
123     ierr = VecNorm(metric, NORM_2, &norm);CHKERRQ(ierr);
124     ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr);
125     errornorm /= norm;
126     if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed (L2 error %f)", errornorm);
127     ierr = VecDestroy(&metricComb);CHKERRQ(ierr);
128 
129     /* Test metric intersection */
130     if (isotropic) {
131       ierr = DMPlexMetricIntersection(dm, 2, metrics, &metricComb);CHKERRQ(ierr);
132       ierr = VecAXPY(metricComb, -1, metric1);CHKERRQ(ierr);
133       ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr);
134       errornorm /= norm;
135       if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed (L2 error %f)", errornorm);
136     }
137     ierr = VecDestroy(&metric2);CHKERRQ(ierr);
138     ierr = VecDestroy(&metricComb);CHKERRQ(ierr);
139     ierr = VecCopy(metric, metric1);CHKERRQ(ierr);
140 
141     /* Test metric SPD enforcement */
142     ierr = DMPlexMetricEnforceSPD(dm, PETSC_TRUE, PETSC_TRUE, metric);CHKERRQ(ierr);
143     if (isotropic) {
144       ierr = VecAXPY(metric1, -1, metric);CHKERRQ(ierr);
145       ierr = VecNorm(metric1, NORM_2, &errornorm);CHKERRQ(ierr);
146       errornorm /= norm;
147       if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed (L2 error %f)", errornorm);
148     }
149     ierr = VecDestroy(&metric1);CHKERRQ(ierr);
150 
151     /* Test metric normalization */
152     ierr = DMPlexMetricNormalize(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1);CHKERRQ(ierr);
153     if (isotropic) {
154       PetscReal target;
155 
156       ierr = DMPlexMetricGetTargetComplexity(dm, &target);CHKERRQ(ierr);
157       scaling = PetscPowReal(target, 2.0/dim);
158       ierr = DMPlexMetricCreateUniform(dm, 0, scaling, &metric2);CHKERRQ(ierr);
159       ierr = VecAXPY(metric2, -1, metric1);CHKERRQ(ierr);
160       ierr = VecNorm(metric2, NORM_2, &errornorm);CHKERRQ(ierr);
161       errornorm /= norm;
162       if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed (L2 error %f)", errornorm);
163     }
164     ierr = VecCopy(metric1, metric);CHKERRQ(ierr);
165     ierr = VecDestroy(&metric2);CHKERRQ(ierr);
166     ierr = VecDestroy(&metric1);CHKERRQ(ierr);
167   }
168 
169   /* Adapt the mesh */
170   ierr = DMAdaptMetric(dm, metric, bdLabel, &dmAdapt);CHKERRQ(ierr);
171   ierr = PetscObjectSetName((PetscObject) dmAdapt, "DM_adapted");CHKERRQ(ierr);
172   ierr = VecDestroy(&metric);CHKERRQ(ierr);
173   ierr = DMViewFromOptions(dmAdapt, NULL, "-adapted_mesh_view");CHKERRQ(ierr);
174 
175   /* Compare DMs */
176   ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
177   ierr = DMView(dmAdapt, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
178 
179   /* Clean up */
180   ierr = DMDestroy(&dmAdapt);CHKERRQ(ierr);
181   ierr = DMDestroy(&dm);CHKERRQ(ierr);
182   ierr = PetscFinalize();
183   return 0;
184 }
185 
186 /*TEST
187 
188   test:
189     suffix: uniform_2d_pragmatic
190     requires: pragmatic
191     args: -dm_plex_metric_target_complexity 100 -dim 2 -dm_adaptor pragmatic -uniform -isotropic
192   test:
193     suffix: uniform_3d_pragmatic
194     requires: pragmatic tetgen
195     args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor pragmatic -uniform -isotropic
196   test:
197     suffix: iso_2d_pragmatic
198     requires: pragmatic
199     args: -dm_plex_metric_target_complexity 100 -dim 2 -dm_adaptor pragmatic -isotropic
200   test:
201     suffix: iso_3d_pragmatic
202     requires: pragmatic tetgen
203     args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor pragmatic -isotropic
204   test:
205     suffix: hessian_2d_pragmatic
206     requires: pragmatic
207     args: -dm_plex_metric_target_complexity 100 -dim 2 -dm_adaptor pragmatic
208   test:
209     suffix: hessian_3d_pragmatic
210     requires: pragmatic tetgen
211     args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor pragmatic
212   test:
213     suffix: uniform_2d_mmg
214     requires: mmg
215     args: -dm_plex_metric_target_complexity 100 -dim 2 -dm_adaptor mmg -uniform -isotropic
216   test:
217     suffix: uniform_3d_mmg
218     requires: mmg tetgen
219     args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor mmg -uniform -isotropic
220   test:
221     suffix: iso_2d_mmg
222     requires: mmg
223     args: -dm_plex_metric_target_complexity 100 -dim 2 -dm_adaptor mmg -isotropic
224   test:
225     suffix: iso_3d_mmg
226     requires: mmg tetgen
227     args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor mmg -isotropic
228   test:
229     suffix: hessian_2d_mmg
230     requires: mmg
231     args: -dm_plex_metric_target_complexity 100 -dim 2 -dm_adaptor mmg
232   test:
233     suffix: hessian_3d_mmg
234     requires: mmg tetgen
235     args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor mmg
236   test:
237     suffix: uniform_3d_parmmg
238     requires: parmmg tetgen
239     nsize: 2
240     args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor parmmg -uniform -isotropic
241   test:
242     suffix: iso_3d_parmmg
243     requires: parmmg tetgen
244     nsize: 2
245     args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor parmmg -isotropic
246   test:
247     suffix: hessian_3d_parmmg
248     requires: parmmg tetgen
249     nsize: 2
250     args: -dm_plex_metric_target_complexity 100 -dim 3 -dm_adaptor parmmg
251 
252 TEST*/
253