xref: /petsc/src/dm/impls/plex/tests/ex19.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown static char help[] = "Tests mesh adaptation with DMPlex and pragmatic.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petsc/private/dmpleximpl.h>
4c4762a1bSJed Brown 
58e3a2eefSMatthew G. Knepley #include <petscsnes.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown typedef struct {
8071b71afSMatthew G. Knepley   PetscInt  Nr;         /* The number of refinement passes */
9c4762a1bSJed Brown   PetscInt  metOpt;     /* Different choices of metric */
10c4762a1bSJed Brown   PetscReal hmax, hmin; /* Max and min sizes prescribed by the metric */
11c4762a1bSJed Brown   PetscBool doL2;       /* Test L2 projection */
12c4762a1bSJed Brown } AppCtx;
13c4762a1bSJed Brown 
14071b71afSMatthew G. Knepley /*
15071b71afSMatthew G. Knepley Classic hyperbolic sensor function for testing multi-scale anisotropic mesh adaptation:
16071b71afSMatthew G. Knepley 
17071b71afSMatthew G. Knepley   f:[-1, 1]x[-1, 1] \to R,
18071b71afSMatthew G. Knepley     f(x, y) = sin(50xy)/100 if |xy| > 2\pi/50 else sin(50xy)
19071b71afSMatthew G. Knepley 
20071b71afSMatthew G. Knepley (mapped to have domain [0,1] x [0,1] in this case).
21071b71afSMatthew G. Knepley */
22071b71afSMatthew G. Knepley static PetscErrorCode sensor(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx)
23071b71afSMatthew G. Knepley {
24071b71afSMatthew G. Knepley   const PetscReal xref = 2.*x[0] - 1.;
25071b71afSMatthew G. Knepley   const PetscReal yref = 2.*x[1] - 1.;
26071b71afSMatthew G. Knepley   const PetscReal xy   = xref*yref;
27071b71afSMatthew G. Knepley 
28071b71afSMatthew G. Knepley   PetscFunctionBeginUser;
29071b71afSMatthew G. Knepley   u[0] = PetscSinReal(50.*xy);
30071b71afSMatthew G. Knepley   if (PetscAbsReal(xy) > 2.*PETSC_PI/50.) u[0] *= 0.01;
31071b71afSMatthew G. Knepley   PetscFunctionReturn(0);
32071b71afSMatthew G. Knepley }
33071b71afSMatthew G. Knepley 
34c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
35c4762a1bSJed Brown {
36c4762a1bSJed Brown   PetscErrorCode ierr;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   PetscFunctionBegin;
39071b71afSMatthew G. Knepley   options->Nr     = 1;
40c4762a1bSJed Brown   options->metOpt = 1;
41c4762a1bSJed Brown   options->hmin   = 0.05;
42c4762a1bSJed Brown   options->hmax   = 0.5;
43c4762a1bSJed Brown   options->doL2   = PETSC_FALSE;
44c4762a1bSJed Brown 
45*9566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(comm, "", "Meshing Adaptation Options", "DMPLEX");PetscCall(ierr);
46*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-Nr", "Numberof refinement passes", "ex19.c", options->Nr, &options->Nr, NULL, 1));
47*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-met", "Different choices of metric", "ex19.c", options->metOpt, &options->metOpt, NULL,0));
48*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-hmax", "Max size prescribed by the metric", "ex19.c", options->hmax, &options->hmax, NULL));
49*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-hmin", "Min size prescribed by the metric", "ex19.c", options->hmin, &options->hmin, NULL));
50*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-do_L2", "Test L2 projection", "ex19.c", options->doL2, &options->doL2, NULL));
51*9566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   PetscFunctionReturn(0);
54c4762a1bSJed Brown }
55c4762a1bSJed Brown 
56071b71afSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
57c4762a1bSJed Brown {
58c4762a1bSJed Brown   PetscFunctionBegin;
59*9566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
60*9566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
61*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
62*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *dm, "DMinit"));
63*9566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-init_dm_view"));
64071b71afSMatthew G. Knepley   PetscFunctionReturn(0);
65c4762a1bSJed Brown }
66071b71afSMatthew G. Knepley 
67071b71afSMatthew G. Knepley static PetscErrorCode ComputeMetricSensor(DM dm, AppCtx *user, Vec *metric)
68c4762a1bSJed Brown {
69071b71afSMatthew G. Knepley   PetscSimplePointFunc funcs[1] = {sensor};
70071b71afSMatthew G. Knepley   DM             dmSensor, dmGrad, dmHess;
71071b71afSMatthew G. Knepley   PetscFE        fe;
72071b71afSMatthew G. Knepley   Vec            f, g, H;
73071b71afSMatthew G. Knepley   PetscBool      simplex;
74071b71afSMatthew G. Knepley   PetscInt       dim;
75c4762a1bSJed Brown 
76071b71afSMatthew G. Knepley   PetscFunctionBegin;
77*9566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
78*9566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
79071b71afSMatthew G. Knepley 
80*9566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dmSensor));
81*9566063dSJacob Faibussowitsch   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, simplex, 1, -1, &fe));
82*9566063dSJacob Faibussowitsch   PetscCall(DMSetField(dmSensor, 0, NULL, (PetscObject) fe));
83*9566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
84*9566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dmSensor));
85*9566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(dmSensor, &f));
86*9566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLocal(dmSensor, 0., funcs, NULL, INSERT_VALUES, f));
87*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(f, NULL, "-sensor_view"));
88071b71afSMatthew G. Knepley 
89071b71afSMatthew G. Knepley   // Recover the gradient of the sensor function
90*9566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dmGrad));
91*9566063dSJacob Faibussowitsch   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, dim, simplex, 1, -1, &fe));
92*9566063dSJacob Faibussowitsch   PetscCall(DMSetField(dmGrad, 0, NULL, (PetscObject) fe));
93*9566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
94*9566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dmGrad));
95*9566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(dmGrad, &g));
96*9566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeGradientClementInterpolant(dmSensor, f, g));
97*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&f));
98*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(g, NULL, "-gradient_view"));
99071b71afSMatthew G. Knepley 
100071b71afSMatthew G. Knepley   // Recover the Hessian of the sensor function
101*9566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dmHess));
102*9566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dmHess, 0, &H));
103*9566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, g, H));
104*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&g));
105*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(H, NULL, "-hessian_view"));
106071b71afSMatthew G. Knepley 
107071b71afSMatthew G. Knepley   // Obtain a metric by Lp normalization
108*9566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricNormalize(dmHess, H, PETSC_TRUE, PETSC_TRUE, metric));
109*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&H));
110*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmHess));
111*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmGrad));
112*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmSensor));
113c4762a1bSJed Brown   PetscFunctionReturn(0);
114c4762a1bSJed Brown }
115c4762a1bSJed Brown 
116c4762a1bSJed Brown static PetscErrorCode ComputeMetric(DM dm, AppCtx *user, Vec *metric)
117c4762a1bSJed Brown {
118071b71afSMatthew G. Knepley   PetscReal          lambda = 1/(user->hmax*user->hmax);
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   PetscFunctionBeginUser;
121071b71afSMatthew G. Knepley   if (user->metOpt == 0) {
1223e32e87aSJoe Wallwork     /* Specify a uniform, isotropic metric */
123*9566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricCreateUniform(dm, 0, lambda, metric));
124071b71afSMatthew G. Knepley   } else if (user->metOpt == 3) {
125*9566063dSJacob Faibussowitsch     PetscCall(ComputeMetricSensor(dm, user, metric));
126071b71afSMatthew G. Knepley   } else {
1273e32e87aSJoe Wallwork     DM                 cdm;
1283e32e87aSJoe Wallwork     Vec                coordinates;
1293e32e87aSJoe Wallwork     const PetscScalar *coords;
1303e32e87aSJoe Wallwork     PetscScalar       *met;
1313e32e87aSJoe Wallwork     PetscReal          h;
132fdc00aefSJoe Wallwork     PetscInt           dim, i, j, vStart, vEnd, v;
1333e32e87aSJoe Wallwork 
134*9566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricCreate(dm, 0, metric));
135*9566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
136*9566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm, &cdm));
137*9566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
138*9566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(coordinates, &coords));
139*9566063dSJacob Faibussowitsch     PetscCall(VecGetArray(*metric, &met));
140*9566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
141071b71afSMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
142071b71afSMatthew G. Knepley       PetscScalar *vcoords;
143c4762a1bSJed Brown       PetscScalar *pmet;
144c4762a1bSJed Brown 
145*9566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRead(cdm, v, coords, &vcoords));
146c4762a1bSJed Brown       switch (user->metOpt) {
147c4762a1bSJed Brown       case 1:
148071b71afSMatthew G. Knepley         h = user->hmax - (user->hmax-user->hmin)*PetscRealPart(vcoords[0]);
149c4762a1bSJed Brown         break;
150c4762a1bSJed Brown       case 2:
151071b71afSMatthew G. Knepley         h = user->hmax*PetscAbsReal(((PetscReal) 1.0)-PetscExpReal(-PetscAbsScalar(vcoords[0]-(PetscReal)0.5))) + user->hmin;
152c4762a1bSJed Brown         break;
153c4762a1bSJed Brown       default:
15498921bdaSJacob Faibussowitsch         SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "metOpt = 0, 1, 2 or 3, cannot be %d", user->metOpt);
155c4762a1bSJed Brown       }
156*9566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &pmet));
157fdc00aefSJoe Wallwork       for (i = 0; i < dim; ++i) {
158fdc00aefSJoe Wallwork         for (j = 0; j < dim; ++j) {
159fdc00aefSJoe Wallwork           if (i == j) {
160fdc00aefSJoe Wallwork             if (i == 0) pmet[i*dim+j] = 1/(h*h);
161fdc00aefSJoe Wallwork             else pmet[i*dim+j] = lambda;
162fdc00aefSJoe Wallwork           } else pmet[i*dim+j] = 0.0;
163fdc00aefSJoe Wallwork         }
164fdc00aefSJoe Wallwork       }
165c4762a1bSJed Brown     }
166*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(*metric, &met));
167*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coordinates, &coords));
1683e32e87aSJoe Wallwork   }
169c4762a1bSJed Brown   PetscFunctionReturn(0);
170c4762a1bSJed Brown }
171c4762a1bSJed Brown 
172c4762a1bSJed Brown static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx)
173c4762a1bSJed Brown {
174c4762a1bSJed Brown   u[0] = x[0] + x[1];
175c4762a1bSJed Brown   return 0;
176c4762a1bSJed Brown }
177c4762a1bSJed Brown 
178c4762a1bSJed Brown static PetscErrorCode TestL2Projection(DM dm, DM dma, AppCtx *user)
179c4762a1bSJed Brown {
180071b71afSMatthew G. Knepley   PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = {linear};
181071b71afSMatthew G. Knepley   DM               dmProj, dmaProj;
182c4762a1bSJed Brown   PetscFE          fe;
183071b71afSMatthew G. Knepley   KSP              ksp;
184071b71afSMatthew G. Knepley   Mat              Interp, mass, mass2;
185071b71afSMatthew G. Knepley   Vec              u, ua, scaling, rhs, uproj;
186c4762a1bSJed Brown   PetscReal        error;
187071b71afSMatthew G. Knepley   PetscBool        simplex;
188c4762a1bSJed Brown   PetscInt         dim;
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   PetscFunctionBeginUser;
191*9566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
192*9566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
193c4762a1bSJed Brown 
194*9566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dmProj));
195*9566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe));
196*9566063dSJacob Faibussowitsch   PetscCall(DMSetField(dmProj, 0, NULL, (PetscObject) fe));
197*9566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
198*9566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dmProj));
199071b71afSMatthew G. Knepley 
200*9566063dSJacob Faibussowitsch   PetscCall(DMClone(dma, &dmaProj));
201*9566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe));
202*9566063dSJacob Faibussowitsch   PetscCall(DMSetField(dmaProj, 0, NULL, (PetscObject) fe));
203*9566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
204*9566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dmaProj));
205071b71afSMatthew G. Knepley 
206*9566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmProj, &u));
207*9566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmaProj, &ua));
208*9566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmaProj, &rhs));
209*9566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmaProj, &uproj));
210071b71afSMatthew G. Knepley 
211071b71afSMatthew G. Knepley   // Interpolate onto original mesh using dual basis
212*9566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dmProj, 0.0, funcs, NULL, INSERT_VALUES, u));
213*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) u, "Original"));
214*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-orig_vec_view"));
215*9566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dmProj, 0.0, funcs, NULL, u, &error));
216*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Original L2 Error: %g\n", (double) error));
217071b71afSMatthew G. Knepley   // Interpolate onto NEW mesh using dual basis
218*9566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dmaProj, 0.0, funcs, NULL, INSERT_VALUES, ua));
219*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) ua, "Adapted"));
220*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(ua, NULL, "-adapt_vec_view"));
221*9566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, ua, &error));
222*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Adapted L2 Error: %g\n", (double) error));
223071b71afSMatthew G. Knepley   // Interpolate between meshes using interpolation matrix
224*9566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(dmProj, dmaProj, &Interp, &scaling));
225*9566063dSJacob Faibussowitsch   PetscCall(MatInterpolate(Interp, u, ua));
226*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Interp));
227*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&scaling));
228*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) ua, "Interpolation"));
229*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(ua, NULL, "-interp_vec_view"));
230*9566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, ua, &error));
231*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Interpolated L2 Error: %g\n", (double) error));
232071b71afSMatthew G. Knepley   // L2 projection
233*9566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(dmaProj, dmaProj, &mass));
234*9566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(mass, NULL, "-mass_mat_view"));
235*9566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
236*9566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, mass, mass));
237*9566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(ksp));
238071b71afSMatthew G. Knepley   //   Compute rhs as M f, could also direclty project the analytic function but we might not have it
239*9566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(dmProj, dmaProj, &mass2));
240*9566063dSJacob Faibussowitsch   PetscCall(MatMult(mass2, u, rhs));
241*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mass2));
242*9566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ksp, rhs, uproj));
243*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) uproj, "L_2 Projection"));
244*9566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(uproj, NULL, "-proj_vec_view"));
245*9566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, uproj, &error));
246*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double) error));
247*9566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ksp));
248*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mass));
249*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmProj, &u));
250*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmaProj, &ua));
251*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmaProj, &rhs));
252*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmaProj, &uproj));
253*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmProj));
254*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmaProj));
255c4762a1bSJed Brown   PetscFunctionReturn(0);
256c4762a1bSJed Brown }
257c4762a1bSJed Brown 
258c4762a1bSJed Brown int main (int argc, char * argv[]) {
259071b71afSMatthew G. Knepley   DM             dm;
260c4762a1bSJed Brown   AppCtx         user;                 /* user-defined work context */
261c4762a1bSJed Brown   MPI_Comm       comm;
262c4762a1bSJed Brown   DM             dma, odm;
263c4762a1bSJed Brown   Vec            metric;
264071b71afSMatthew G. Knepley   PetscInt       r;
265c4762a1bSJed Brown 
266*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
267c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
268*9566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(comm, &user));
269*9566063dSJacob Faibussowitsch   PetscCall(CreateMesh(comm, &dm));
270c4762a1bSJed Brown 
271071b71afSMatthew G. Knepley   odm  = dm;
272*9566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeOverlap(odm, 1, NULL, &dm));
273071b71afSMatthew G. Knepley   if (!dm) {dm = odm;}
274*9566063dSJacob Faibussowitsch   else     PetscCall(DMDestroy(&odm));
275071b71afSMatthew G. Knepley 
276071b71afSMatthew G. Knepley   for (r = 0; r < user.Nr; ++r) {
277071b71afSMatthew G. Knepley     DMLabel label;
278071b71afSMatthew G. Knepley 
279*9566063dSJacob Faibussowitsch     PetscCall(ComputeMetric(dm, &user, &metric));
280*9566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, "marker", &label));
281*9566063dSJacob Faibussowitsch     PetscCall(DMAdaptMetric(dm, metric, label, NULL, &dma));
282*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&metric));
283*9566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) dma, "DMadapt"));
284*9566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject) dma, "adapt_"));
285*9566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dma, NULL, "-dm_view"));
286*9566063dSJacob Faibussowitsch     if (user.doL2) PetscCall(TestL2Projection(dm, dma, &user));
287*9566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dm));
288071b71afSMatthew G. Knepley     dm   = dma;
289071b71afSMatthew G. Knepley   }
290*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) dm, "final_"));
291*9566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
292*9566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
293*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
294b122ec5aSJacob Faibussowitsch   return 0;
295c4762a1bSJed Brown }
296c4762a1bSJed Brown 
297c4762a1bSJed Brown /*TEST
298c4762a1bSJed Brown 
2998e3a2eefSMatthew G. Knepley   build:
3008e3a2eefSMatthew G. Knepley     requires: pragmatic
3018e3a2eefSMatthew G. Knepley 
302071b71afSMatthew G. Knepley   testset:
3039e3182e6SJoe Wallwork     args: -dm_plex_box_faces 4,4,4 -dm_adaptor pragmatic -met 2 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic
304071b71afSMatthew G. Knepley 
305c4762a1bSJed Brown     test:
306c4762a1bSJed Brown       suffix: 0
307071b71afSMatthew G. Knepley       args: -dm_plex_separate_marker 0
308c4762a1bSJed Brown     test:
309c4762a1bSJed Brown       suffix: 1
310071b71afSMatthew G. Knepley       args: -dm_plex_separate_marker 1
311c4762a1bSJed Brown     test:
312c4762a1bSJed Brown       suffix: 2
313071b71afSMatthew G. Knepley       args: -dm_plex_dim 3
314c4762a1bSJed Brown     test:
315c4762a1bSJed Brown       suffix: 3
316071b71afSMatthew G. Knepley       args: -dm_plex_dim 3
317071b71afSMatthew G. Knepley 
318071b71afSMatthew G. Knepley   # Pragmatic hangs for simple partitioner
319071b71afSMatthew G. Knepley   testset:
320071b71afSMatthew G. Knepley     requires: parmetis
321e600fa54SMatthew G. Knepley     args: -dm_plex_box_faces 2,2 -dm_adaptor pragmatic -petscpartitioner_type parmetis -met 2 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic
322071b71afSMatthew G. Knepley 
323c4762a1bSJed Brown     test:
324c4762a1bSJed Brown       suffix: 4
325c4762a1bSJed Brown       nsize: 2
326c4762a1bSJed Brown     test:
327c4762a1bSJed Brown       suffix: 5
328c4762a1bSJed Brown       nsize: 4
329071b71afSMatthew G. Knepley 
330c4762a1bSJed Brown   test:
3319880b877SBarry Smith     requires: parmetis
332c4762a1bSJed Brown     suffix: 6
333c4762a1bSJed Brown     nsize: 2
334e600fa54SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 9,9,9 -dm_adaptor pragmatic -petscpartitioner_type parmetis \
3359e3182e6SJoe Wallwork           -met 0 -hmin 0.01 -hmax 0.03 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic
336c4762a1bSJed Brown   test:
3379880b877SBarry Smith     requires: parmetis
338c4762a1bSJed Brown     suffix: 7
3393e32e87aSJoe Wallwork     nsize: 2
340e600fa54SMatthew G. Knepley     args: -dm_plex_box_faces 19,19 -dm_adaptor pragmatic -petscpartitioner_type parmetis \
3419e3182e6SJoe Wallwork           -met 2 -hmax 0.5 -hmin 0.001 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic
342c4762a1bSJed Brown   test:
343c4762a1bSJed Brown     suffix: proj_0
344cbddfcd8SJoe Wallwork     args: -dm_plex_box_faces 2,2 -dm_plex_hash_location -dm_adaptor pragmatic -init_dm_view -adapt_dm_view -do_L2 \
3459e3182e6SJoe Wallwork           -petscspace_degree 1 -petscfe_default_quadrature_order 1 -pc_type lu -dm_adaptor pragmatic
346c4762a1bSJed Brown   test:
347c4762a1bSJed Brown     suffix: proj_1
348cbddfcd8SJoe Wallwork     args: -dm_plex_box_faces 4,4 -dm_plex_hash_location -dm_adaptor pragmatic -init_dm_view -adapt_dm_view -do_L2 \
3499e3182e6SJoe Wallwork           -petscspace_degree 1 -petscfe_default_quadrature_order 1 -pc_type lu -dm_adaptor pragmatic
350071b71afSMatthew G. Knepley 
351071b71afSMatthew G. Knepley   test:
352071b71afSMatthew G. Knepley     suffix: sensor
353cbddfcd8SJoe Wallwork     args: -dm_plex_box_faces 9,9 -met 3 -dm_adaptor pragmatic -init_dm_view -adapt_dm_view \
354071b71afSMatthew G. Knepley           -dm_plex_metric_h_min 1.e-10 -dm_plex_metric_h_max 1.0e-01 -dm_plex_metric_a_max 1.0e+05 -dm_plex_metric_p 1.0 \
3559e3182e6SJoe Wallwork             -dm_plex_metric_target_complexity 10000.0 -dm_adaptor pragmatic
356c4762a1bSJed Brown 
357c4762a1bSJed Brown TEST*/
358