xref: /petsc/src/dm/impls/plex/tests/ex60.c (revision a4af0ceea8a251db97ee0dc5c0d52d4adf50264a)
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 = PetscOptionsInt("-num_edges", "Number of edges on each boundary of the initial mesh", "ex60.c", numEdges, &numEdges, NULL);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 = DMAddField(dmIndi, 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       PetscInt         pStart, pEnd, p;
80       PetscSection     sec, secCoord;
81       Vec              gradient;
82 
83       /* Project the parabola into P1 space */
84       ierr = DMProjectFunctionLocal(dmIndi, 0.0, funcs, NULL, INSERT_ALL_VALUES, indicator);CHKERRQ(ierr);
85 
86       /* Approximate the gradient */
87       ierr = DMClone(dmIndi, &dmGrad);CHKERRQ(ierr);
88       ierr = DMGetCoordinateSection(dmGrad, &secCoord);CHKERRQ(ierr);
89       ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dmGrad), &sec);CHKERRQ(ierr);
90       ierr = PetscSectionSetNumFields(sec, 1);CHKERRQ(ierr);
91       ierr = PetscSectionSetFieldComponents(sec, 0, dim);CHKERRQ(ierr);
92       ierr = PetscSectionGetChart(secCoord, &pStart, &pEnd);CHKERRQ(ierr);
93       ierr = PetscSectionSetChart(sec, pStart, pEnd);CHKERRQ(ierr);
94       for (p = pStart; p < pEnd; ++p) {
95         ierr = PetscSectionSetDof(sec, p, dim);CHKERRQ(ierr);
96         ierr = PetscSectionSetFieldDof(sec, p, 0, dim);CHKERRQ(ierr);
97       }
98       ierr = PetscSectionSetUp(sec);CHKERRQ(ierr);
99       ierr = DMSetLocalSection(dmGrad, sec);CHKERRQ(ierr);
100       ierr = PetscSectionDestroy(&sec);CHKERRQ(ierr);
101       ierr = PetscFECreateLagrange(comm, dim, dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
102       ierr = DMClearFields(dmGrad);CHKERRQ(ierr);
103       ierr = DMAddField(dmGrad, NULL, (PetscObject)fe);CHKERRQ(ierr);
104       ierr = DMCreateDS(dmGrad);CHKERRQ(ierr);
105       ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
106       ierr = DMCreateLocalVector(dmGrad, &gradient);CHKERRQ(ierr);
107       ierr = DMPlexComputeGradientClementInterpolant(dmIndi, indicator, gradient);CHKERRQ(ierr);
108       ierr = VecViewFromOptions(gradient, NULL, "-adapt_gradient_view");CHKERRQ(ierr);
109 
110       /* Approximate the Hessian */
111       ierr = DMGetCoordinateSection(dm, &secCoord);CHKERRQ(ierr);
112       ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &sec);CHKERRQ(ierr);
113       ierr = PetscSectionSetNumFields(sec, 1);CHKERRQ(ierr);
114       ierr = PetscSectionSetFieldComponents(sec, 0, dim*dim);CHKERRQ(ierr);
115       ierr = PetscSectionGetChart(secCoord, &pStart, &pEnd);CHKERRQ(ierr);
116       ierr = PetscSectionSetChart(sec, pStart, pEnd);CHKERRQ(ierr);
117       for (p = pStart; p < pEnd; ++p) {
118         ierr = PetscSectionSetDof(sec, p, dim*dim);CHKERRQ(ierr);
119         ierr = PetscSectionSetFieldDof(sec, p, 0, dim*dim);CHKERRQ(ierr);
120       }
121       ierr = PetscSectionSetUp(sec);CHKERRQ(ierr);
122       ierr = DMSetLocalSection(dm, sec);CHKERRQ(ierr);
123       ierr = PetscSectionDestroy(&sec);CHKERRQ(ierr);
124       ierr = PetscFECreateLagrange(comm, dim, dim*dim, PETSC_TRUE, 1, PETSC_DETERMINE, &fe);CHKERRQ(ierr);
125       ierr = DMClearFields(dm);CHKERRQ(ierr);
126       ierr = DMAddField(dm, NULL, (PetscObject)fe);CHKERRQ(ierr);
127       ierr = DMCreateDS(dm);CHKERRQ(ierr);
128       ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
129       ierr = DMCreateLocalVector(dm, &metric);CHKERRQ(ierr);
130       ierr = DMPlexComputeGradientClementInterpolant(dmGrad, gradient, metric);CHKERRQ(ierr);
131       ierr = VecViewFromOptions(metric, NULL, "-adapt_hessian_view");CHKERRQ(ierr);
132       ierr = VecDestroy(&gradient);CHKERRQ(ierr);
133       ierr = DMDestroy(&dmGrad);CHKERRQ(ierr);
134     }
135     ierr = VecDestroy(&indicator);CHKERRQ(ierr);
136     ierr = DMDestroy(&dmIndi);CHKERRQ(ierr);
137   }
138 
139   /* Test metric routines */
140   {
141     PetscReal   errornorm, norm, tol = 1.0e-10, weights[2] = {0.8, 0.2};
142     Vec         metric1, metric2, metricComb;
143     Vec         metrics[2];
144 
145     ierr = VecDuplicate(metric, &metric1);CHKERRQ(ierr);
146     ierr = VecSet(metric1, 0);CHKERRQ(ierr);
147     ierr = VecAXPY(metric1, 0.625, metric);CHKERRQ(ierr);
148     ierr = VecDuplicate(metric, &metric2);CHKERRQ(ierr);
149     ierr = VecSet(metric2, 0);CHKERRQ(ierr);
150     ierr = VecAXPY(metric2, 2.5, metric);CHKERRQ(ierr);
151     metrics[0] = metric1;
152     metrics[1] = metric2;
153 
154     /* Test metric average */
155     ierr = DMPlexMetricAverage(dm, 2, weights, metrics, &metricComb);CHKERRQ(ierr);
156     ierr = VecAXPY(metricComb, -1, metric);CHKERRQ(ierr);
157     ierr = VecNorm(metric, NORM_2, &norm);CHKERRQ(ierr);
158     ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr);
159     errornorm /= norm;
160     if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric average test failed (L2 error %f)", errornorm);
161     ierr = VecDestroy(&metricComb);CHKERRQ(ierr);
162 
163     /* Test metric intersection */
164     if (isotropic) {
165       ierr = DMPlexMetricIntersection(dm, 2, metrics, &metricComb);CHKERRQ(ierr);
166       ierr = VecAXPY(metricComb, -1, metric1);CHKERRQ(ierr);
167       ierr = VecNorm(metricComb, NORM_2, &errornorm);CHKERRQ(ierr);
168       errornorm /= norm;
169       if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric intersection test failed (L2 error %f)", errornorm);
170     }
171     ierr = VecDestroy(&metric2);CHKERRQ(ierr);
172     ierr = VecDestroy(&metricComb);CHKERRQ(ierr);
173     ierr = VecCopy(metric, metric1);CHKERRQ(ierr);
174 
175     /* Test metric SPD enforcement */
176     ierr = DMPlexMetricEnforceSPD(dm, PETSC_TRUE, PETSC_TRUE, metric);CHKERRQ(ierr);
177     if (isotropic) {
178       ierr = VecAXPY(metric1, -1, metric);CHKERRQ(ierr);
179       ierr = VecNorm(metric1, NORM_2, &errornorm);CHKERRQ(ierr);
180       errornorm /= norm;
181       if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric SPD enforcement test failed (L2 error %f)", errornorm);
182     }
183     ierr = VecDestroy(&metric1);CHKERRQ(ierr);
184 
185     /* Test metric normalization */
186     ierr = DMPlexMetricNormalize(dm, metric, PETSC_TRUE, PETSC_TRUE, &metric1);CHKERRQ(ierr);
187     if (isotropic) {
188       PetscReal target;
189 
190       ierr = DMPlexMetricGetTargetComplexity(dm, &target);CHKERRQ(ierr);
191       scaling = PetscPowReal(target, 2.0/dim);
192       ierr = DMPlexMetricCreateUniform(dm, 0, scaling, &metric2);CHKERRQ(ierr);
193       ierr = VecAXPY(metric2, -1, metric1);CHKERRQ(ierr);
194       ierr = VecNorm(metric2, NORM_2, &errornorm);CHKERRQ(ierr);
195       errornorm /= norm;
196       if (errornorm > tol) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Metric normalization test failed (L2 error %f)", errornorm);
197     }
198     ierr = VecCopy(metric1, metric);CHKERRQ(ierr);
199     ierr = VecDestroy(&metric2);CHKERRQ(ierr);
200     ierr = VecDestroy(&metric1);CHKERRQ(ierr);
201   }
202 
203   /* Adapt the mesh */
204   ierr = DMAdaptMetric(dm, metric, bdLabel, &dmAdapt);CHKERRQ(ierr);
205   ierr = PetscObjectSetName((PetscObject) dmAdapt, "DM_adapted");CHKERRQ(ierr);
206   ierr = VecDestroy(&metric);CHKERRQ(ierr);
207   ierr = DMViewFromOptions(dmAdapt, NULL, "-adapted_mesh_view");CHKERRQ(ierr);
208 
209   /* Compare DMs */
210   ierr = DMView(dm, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
211   ierr = DMView(dmAdapt, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
212 
213   /* Clean up */
214   ierr = DMDestroy(&dmAdapt);CHKERRQ(ierr);
215   ierr = DMDestroy(&dm);CHKERRQ(ierr);
216   ierr = PetscFinalize();
217   return 0;
218 }
219 
220 /*TEST
221 
222   build:
223     requires: pragmatic
224 
225   test:
226     suffix: uniform_2d
227     args: -dm_plex_metric_target_complexity 100 -dim 2 -uniform -isotropic
228   test:
229     suffix: uniform_3d
230     args: -dm_plex_metric_target_complexity 100 -dim 3 -uniform -isotropic
231   test:
232     suffix: iso_2d
233     args: -dm_plex_metric_target_complexity 100 -dim 2 -isotropic
234   test:
235     suffix: iso_3d
236     args: -dm_plex_metric_target_complexity 100 -dim 3 -isotropic
237   test:
238     suffix: hessian_2d
239     args: -dm_plex_metric_target_complexity 100 -dim 2
240   test:
241     suffix: hessian_3d
242     args: -dm_plex_metric_target_complexity 100 -dim 3
243 
244 TEST*/
245