xref: /petsc/src/dm/impls/plex/tests/ex19.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 */
22*9371c9d4SSatish Balay static PetscErrorCode sensor(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) {
23071b71afSMatthew G. Knepley   const PetscReal xref = 2. * x[0] - 1.;
24071b71afSMatthew G. Knepley   const PetscReal yref = 2. * x[1] - 1.;
25071b71afSMatthew G. Knepley   const PetscReal xy   = xref * yref;
26071b71afSMatthew G. Knepley 
27071b71afSMatthew G. Knepley   PetscFunctionBeginUser;
28071b71afSMatthew G. Knepley   u[0] = PetscSinReal(50. * xy);
29071b71afSMatthew G. Knepley   if (PetscAbsReal(xy) > 2. * PETSC_PI / 50.) u[0] *= 0.01;
30071b71afSMatthew G. Knepley   PetscFunctionReturn(0);
31071b71afSMatthew G. Knepley }
32071b71afSMatthew G. Knepley 
33*9371c9d4SSatish Balay static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
34c4762a1bSJed Brown   PetscFunctionBegin;
35071b71afSMatthew G. Knepley   options->Nr     = 1;
36c4762a1bSJed Brown   options->metOpt = 1;
37c4762a1bSJed Brown   options->hmin   = 0.05;
38c4762a1bSJed Brown   options->hmax   = 0.5;
39c4762a1bSJed Brown   options->doL2   = PETSC_FALSE;
40c4762a1bSJed Brown 
41d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Adaptation Options", "DMPLEX");
429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-Nr", "Numberof refinement passes", "ex19.c", options->Nr, &options->Nr, NULL, 1));
439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-met", "Different choices of metric", "ex19.c", options->metOpt, &options->metOpt, NULL, 0));
449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-hmax", "Max size prescribed by the metric", "ex19.c", options->hmax, &options->hmax, NULL));
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-hmin", "Min size prescribed by the metric", "ex19.c", options->hmin, &options->hmin, NULL));
469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-do_L2", "Test L2 projection", "ex19.c", options->doL2, &options->doL2, NULL));
47d0609cedSBarry Smith   PetscOptionsEnd();
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   PetscFunctionReturn(0);
50c4762a1bSJed Brown }
51c4762a1bSJed Brown 
52*9371c9d4SSatish Balay static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) {
53c4762a1bSJed Brown   PetscFunctionBegin;
549566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
559566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
569566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
579566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*dm, "DMinit"));
589566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-init_dm_view"));
59071b71afSMatthew G. Knepley   PetscFunctionReturn(0);
60c4762a1bSJed Brown }
61071b71afSMatthew G. Knepley 
62*9371c9d4SSatish Balay static PetscErrorCode ComputeMetricSensor(DM dm, AppCtx *user, Vec *metric) {
63071b71afSMatthew G. Knepley   PetscSimplePointFunc funcs[1] = {sensor};
64edaa2528SJoe Wallwork   DM                   dmSensor, dmGrad, dmHess, dmDet;
65071b71afSMatthew G. Knepley   PetscFE              fe;
66edaa2528SJoe Wallwork   Vec                  f, g, H, determinant;
67071b71afSMatthew G. Knepley   PetscBool            simplex;
68071b71afSMatthew G. Knepley   PetscInt             dim;
69c4762a1bSJed Brown 
70071b71afSMatthew G. Knepley   PetscFunctionBegin;
719566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
729566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
73071b71afSMatthew G. Knepley 
749566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dmSensor));
759566063dSJacob Faibussowitsch   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, simplex, 1, -1, &fe));
769566063dSJacob Faibussowitsch   PetscCall(DMSetField(dmSensor, 0, NULL, (PetscObject)fe));
779566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
789566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dmSensor));
799566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(dmSensor, &f));
809566063dSJacob Faibussowitsch   PetscCall(DMProjectFunctionLocal(dmSensor, 0., funcs, NULL, INSERT_VALUES, f));
819566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(f, NULL, "-sensor_view"));
82071b71afSMatthew G. Knepley 
83071b71afSMatthew G. Knepley   // Recover the gradient of the sensor function
849566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dmGrad));
859566063dSJacob Faibussowitsch   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, dim, simplex, 1, -1, &fe));
869566063dSJacob Faibussowitsch   PetscCall(DMSetField(dmGrad, 0, NULL, (PetscObject)fe));
879566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
889566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dmGrad));
899566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(dmGrad, &g));
909566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeGradientClementInterpolant(dmSensor, f, g));
919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&f));
929566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(g, NULL, "-gradient_view"));
93071b71afSMatthew G. Knepley 
94071b71afSMatthew G. Knepley   // Recover the Hessian of the sensor function
959566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dmHess));
969566063dSJacob Faibussowitsch   PetscCall(DMPlexMetricCreate(dmHess, 0, &H));
979566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, g, H));
989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&g));
999566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(H, NULL, "-hessian_view"));
100071b71afSMatthew G. Knepley 
101071b71afSMatthew G. Knepley   // Obtain a metric by Lp normalization
102edaa2528SJoe Wallwork   PetscCall(DMPlexMetricCreate(dm, 0, metric));
103edaa2528SJoe Wallwork   PetscCall(DMPlexMetricDeterminantCreate(dm, 0, &determinant, &dmDet));
104edaa2528SJoe Wallwork   PetscCall(DMPlexMetricNormalize(dmHess, H, PETSC_TRUE, PETSC_TRUE, *metric, determinant));
105edaa2528SJoe Wallwork   PetscCall(VecDestroy(&determinant));
106edaa2528SJoe Wallwork   PetscCall(DMDestroy(&dmDet));
1079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&H));
1089566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmHess));
1099566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmGrad));
1109566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmSensor));
111c4762a1bSJed Brown   PetscFunctionReturn(0);
112c4762a1bSJed Brown }
113c4762a1bSJed Brown 
114*9371c9d4SSatish Balay static PetscErrorCode ComputeMetric(DM dm, AppCtx *user, Vec *metric) {
115071b71afSMatthew G. Knepley   PetscReal lambda = 1 / (user->hmax * user->hmax);
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   PetscFunctionBeginUser;
118071b71afSMatthew G. Knepley   if (user->metOpt == 0) {
1193e32e87aSJoe Wallwork     /* Specify a uniform, isotropic metric */
1209566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricCreateUniform(dm, 0, lambda, metric));
121071b71afSMatthew G. Knepley   } else if (user->metOpt == 3) {
1229566063dSJacob Faibussowitsch     PetscCall(ComputeMetricSensor(dm, user, metric));
123071b71afSMatthew G. Knepley   } else {
1243e32e87aSJoe Wallwork     DM                 cdm;
1253e32e87aSJoe Wallwork     Vec                coordinates;
1263e32e87aSJoe Wallwork     const PetscScalar *coords;
1273e32e87aSJoe Wallwork     PetscScalar       *met;
1283e32e87aSJoe Wallwork     PetscReal          h;
129fdc00aefSJoe Wallwork     PetscInt           dim, i, j, vStart, vEnd, v;
1303e32e87aSJoe Wallwork 
1319566063dSJacob Faibussowitsch     PetscCall(DMPlexMetricCreate(dm, 0, metric));
1329566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
1339566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm, &cdm));
1349566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
1359566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(coordinates, &coords));
1369566063dSJacob Faibussowitsch     PetscCall(VecGetArray(*metric, &met));
1379566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
138071b71afSMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
139071b71afSMatthew G. Knepley       PetscScalar *vcoords;
140c4762a1bSJed Brown       PetscScalar *pmet;
141c4762a1bSJed Brown 
1429566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRead(cdm, v, coords, &vcoords));
143c4762a1bSJed Brown       switch (user->metOpt) {
144*9371c9d4SSatish Balay       case 1: h = user->hmax - (user->hmax - user->hmin) * PetscRealPart(vcoords[0]); break;
145*9371c9d4SSatish Balay       case 2: h = user->hmax * PetscAbsReal(((PetscReal)1.0) - PetscExpReal(-PetscAbsScalar(vcoords[0] - (PetscReal)0.5))) + user->hmin; break;
146*9371c9d4SSatish Balay       default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "metOpt = 0, 1, 2 or 3, cannot be %d", user->metOpt);
147c4762a1bSJed Brown       }
1489566063dSJacob Faibussowitsch       PetscCall(DMPlexPointLocalRef(dm, v, met, &pmet));
149fdc00aefSJoe Wallwork       for (i = 0; i < dim; ++i) {
150fdc00aefSJoe Wallwork         for (j = 0; j < dim; ++j) {
151fdc00aefSJoe Wallwork           if (i == j) {
152fdc00aefSJoe Wallwork             if (i == 0) pmet[i * dim + j] = 1 / (h * h);
153fdc00aefSJoe Wallwork             else pmet[i * dim + j] = lambda;
154fdc00aefSJoe Wallwork           } else pmet[i * dim + j] = 0.0;
155fdc00aefSJoe Wallwork         }
156fdc00aefSJoe Wallwork       }
157c4762a1bSJed Brown     }
1589566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(*metric, &met));
1599566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(coordinates, &coords));
1603e32e87aSJoe Wallwork   }
161c4762a1bSJed Brown   PetscFunctionReturn(0);
162c4762a1bSJed Brown }
163c4762a1bSJed Brown 
164*9371c9d4SSatish Balay static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) {
165c4762a1bSJed Brown   u[0] = x[0] + x[1];
166c4762a1bSJed Brown   return 0;
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
169*9371c9d4SSatish Balay static PetscErrorCode TestL2Projection(DM dm, DM dma, AppCtx *user) {
170071b71afSMatthew G. Knepley   PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {linear};
171071b71afSMatthew G. Knepley   DM        dmProj, dmaProj;
172c4762a1bSJed Brown   PetscFE   fe;
173071b71afSMatthew G. Knepley   KSP       ksp;
174071b71afSMatthew G. Knepley   Mat       Interp, mass, mass2;
175071b71afSMatthew G. Knepley   Vec       u, ua, scaling, rhs, uproj;
176c4762a1bSJed Brown   PetscReal error;
177071b71afSMatthew G. Knepley   PetscBool simplex;
178c4762a1bSJed Brown   PetscInt  dim;
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   PetscFunctionBeginUser;
1819566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1829566063dSJacob Faibussowitsch   PetscCall(DMPlexIsSimplex(dm, &simplex));
183c4762a1bSJed Brown 
1849566063dSJacob Faibussowitsch   PetscCall(DMClone(dm, &dmProj));
1859566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe));
1869566063dSJacob Faibussowitsch   PetscCall(DMSetField(dmProj, 0, NULL, (PetscObject)fe));
1879566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
1889566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dmProj));
189071b71afSMatthew G. Knepley 
1909566063dSJacob Faibussowitsch   PetscCall(DMClone(dma, &dmaProj));
1919566063dSJacob Faibussowitsch   PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe));
1929566063dSJacob Faibussowitsch   PetscCall(DMSetField(dmaProj, 0, NULL, (PetscObject)fe));
1939566063dSJacob Faibussowitsch   PetscCall(PetscFEDestroy(&fe));
1949566063dSJacob Faibussowitsch   PetscCall(DMCreateDS(dmaProj));
195071b71afSMatthew G. Knepley 
1969566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmProj, &u));
1979566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmaProj, &ua));
1989566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmaProj, &rhs));
1999566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalVector(dmaProj, &uproj));
200071b71afSMatthew G. Knepley 
201071b71afSMatthew G. Knepley   // Interpolate onto original mesh using dual basis
2029566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dmProj, 0.0, funcs, NULL, INSERT_VALUES, u));
2039566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)u, "Original"));
2049566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(u, NULL, "-orig_vec_view"));
2059566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dmProj, 0.0, funcs, NULL, u, &error));
2069566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Original L2 Error: %g\n", (double)error));
207071b71afSMatthew G. Knepley   // Interpolate onto NEW mesh using dual basis
2089566063dSJacob Faibussowitsch   PetscCall(DMProjectFunction(dmaProj, 0.0, funcs, NULL, INSERT_VALUES, ua));
2099566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)ua, "Adapted"));
2109566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(ua, NULL, "-adapt_vec_view"));
2119566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, ua, &error));
2129566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Adapted L2 Error: %g\n", (double)error));
213071b71afSMatthew G. Knepley   // Interpolate between meshes using interpolation matrix
2149566063dSJacob Faibussowitsch   PetscCall(DMCreateInterpolation(dmProj, dmaProj, &Interp, &scaling));
2159566063dSJacob Faibussowitsch   PetscCall(MatInterpolate(Interp, u, ua));
2169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Interp));
2179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&scaling));
2189566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)ua, "Interpolation"));
2199566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(ua, NULL, "-interp_vec_view"));
2209566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, ua, &error));
2219566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Interpolated L2 Error: %g\n", (double)error));
222071b71afSMatthew G. Knepley   // L2 projection
2239566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(dmaProj, dmaProj, &mass));
2249566063dSJacob Faibussowitsch   PetscCall(MatViewFromOptions(mass, NULL, "-mass_mat_view"));
2259566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
2269566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(ksp, mass, mass));
2279566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(ksp));
228071b71afSMatthew G. Knepley   //   Compute rhs as M f, could also direclty project the analytic function but we might not have it
2299566063dSJacob Faibussowitsch   PetscCall(DMCreateMassMatrix(dmProj, dmaProj, &mass2));
2309566063dSJacob Faibussowitsch   PetscCall(MatMult(mass2, u, rhs));
2319566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mass2));
2329566063dSJacob Faibussowitsch   PetscCall(KSPSolve(ksp, rhs, uproj));
2339566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)uproj, "L_2 Projection"));
2349566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(uproj, NULL, "-proj_vec_view"));
2359566063dSJacob Faibussowitsch   PetscCall(DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, uproj, &error));
2369566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double)error));
2379566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&ksp));
2389566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mass));
2399566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmProj, &u));
2409566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmaProj, &ua));
2419566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmaProj, &rhs));
2429566063dSJacob Faibussowitsch   PetscCall(DMRestoreGlobalVector(dmaProj, &uproj));
2439566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmProj));
2449566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmaProj));
245c4762a1bSJed Brown   PetscFunctionReturn(0);
246c4762a1bSJed Brown }
247c4762a1bSJed Brown 
248c4762a1bSJed Brown int main(int argc, char *argv[]) {
249071b71afSMatthew G. Knepley   DM       dm;
250c4762a1bSJed Brown   AppCtx   user; /* user-defined work context */
251c4762a1bSJed Brown   MPI_Comm comm;
252c4762a1bSJed Brown   DM       dma, odm;
253c4762a1bSJed Brown   Vec      metric;
254071b71afSMatthew G. Knepley   PetscInt r;
255c4762a1bSJed Brown 
256327415f7SBarry Smith   PetscFunctionBeginUser;
2579566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
258c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
2599566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(comm, &user));
2609566063dSJacob Faibussowitsch   PetscCall(CreateMesh(comm, &dm));
261c4762a1bSJed Brown 
262071b71afSMatthew G. Knepley   odm = dm;
2639566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeOverlap(odm, 1, NULL, &dm));
264*9371c9d4SSatish Balay   if (!dm) {
265*9371c9d4SSatish Balay     dm = odm;
266*9371c9d4SSatish Balay   } else PetscCall(DMDestroy(&odm));
267071b71afSMatthew G. Knepley 
268071b71afSMatthew G. Knepley   for (r = 0; r < user.Nr; ++r) {
269071b71afSMatthew G. Knepley     DMLabel label;
270071b71afSMatthew G. Knepley 
2719566063dSJacob Faibussowitsch     PetscCall(ComputeMetric(dm, &user, &metric));
2729566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, "marker", &label));
2739566063dSJacob Faibussowitsch     PetscCall(DMAdaptMetric(dm, metric, label, NULL, &dma));
2749566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&metric));
2759566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)dma, "DMadapt"));
2769566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dma, "adapt_"));
2779566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dma, NULL, "-dm_view"));
2789566063dSJacob Faibussowitsch     if (user.doL2) PetscCall(TestL2Projection(dm, dma, &user));
2799566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dm));
280071b71afSMatthew G. Knepley     dm = dma;
281071b71afSMatthew G. Knepley   }
2829566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, "final_"));
2839566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
2849566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
2859566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
286b122ec5aSJacob Faibussowitsch   return 0;
287c4762a1bSJed Brown }
288c4762a1bSJed Brown 
289c4762a1bSJed Brown /*TEST
290c4762a1bSJed Brown 
2918e3a2eefSMatthew G. Knepley   build:
2928e3a2eefSMatthew G. Knepley     requires: pragmatic
2938e3a2eefSMatthew G. Knepley 
294071b71afSMatthew G. Knepley   testset:
2959e3182e6SJoe Wallwork     args: -dm_plex_box_faces 4,4,4 -dm_adaptor pragmatic -met 2 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic
296071b71afSMatthew G. Knepley 
297c4762a1bSJed Brown     test:
298edaa2528SJoe Wallwork       suffix: 2d
299071b71afSMatthew G. Knepley       args: -dm_plex_separate_marker 0
300c4762a1bSJed Brown     test:
301edaa2528SJoe Wallwork       suffix: 2d_sep
302071b71afSMatthew G. Knepley       args: -dm_plex_separate_marker 1
303c4762a1bSJed Brown     test:
304edaa2528SJoe Wallwork       suffix: 3d
305071b71afSMatthew G. Knepley       args: -dm_plex_dim 3
306071b71afSMatthew G. Knepley 
307071b71afSMatthew G. Knepley   # Pragmatic hangs for simple partitioner
308071b71afSMatthew G. Knepley   testset:
309071b71afSMatthew G. Knepley     requires: parmetis
310edaa2528SJoe Wallwork     args: -dm_plex_box_faces 2,2 -petscpartitioner_type parmetis -met 2 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic
311071b71afSMatthew G. Knepley 
312c4762a1bSJed Brown     test:
313edaa2528SJoe Wallwork       suffix: 2d_parmetis_np2
314c4762a1bSJed Brown       nsize: 2
315c4762a1bSJed Brown     test:
316edaa2528SJoe Wallwork       suffix: 2d_parmetis_np4
317c4762a1bSJed Brown       nsize: 4
318071b71afSMatthew G. Knepley 
319c4762a1bSJed Brown   test:
3209880b877SBarry Smith     requires: parmetis
321edaa2528SJoe Wallwork     suffix: 3d_parmetis_met0
322c4762a1bSJed Brown     nsize: 2
323e600fa54SMatthew G. Knepley     args: -dm_plex_dim 3 -dm_plex_box_faces 9,9,9 -dm_adaptor pragmatic -petscpartitioner_type parmetis \
3249e3182e6SJoe Wallwork           -met 0 -hmin 0.01 -hmax 0.03 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic
325c4762a1bSJed Brown   test:
3269880b877SBarry Smith     requires: parmetis
327edaa2528SJoe Wallwork     suffix: 3d_parmetis_met2
3283e32e87aSJoe Wallwork     nsize: 2
329e600fa54SMatthew G. Knepley     args: -dm_plex_box_faces 19,19 -dm_adaptor pragmatic -petscpartitioner_type parmetis \
3309e3182e6SJoe Wallwork           -met 2 -hmax 0.5 -hmin 0.001 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic
331c4762a1bSJed Brown   test:
332edaa2528SJoe Wallwork     suffix: proj2
333cbddfcd8SJoe Wallwork     args: -dm_plex_box_faces 2,2 -dm_plex_hash_location -dm_adaptor pragmatic -init_dm_view -adapt_dm_view -do_L2 \
3349e3182e6SJoe Wallwork           -petscspace_degree 1 -petscfe_default_quadrature_order 1 -pc_type lu -dm_adaptor pragmatic
335c4762a1bSJed Brown   test:
336edaa2528SJoe Wallwork     suffix: proj4
337cbddfcd8SJoe Wallwork     args: -dm_plex_box_faces 4,4 -dm_plex_hash_location -dm_adaptor pragmatic -init_dm_view -adapt_dm_view -do_L2 \
3389e3182e6SJoe Wallwork           -petscspace_degree 1 -petscfe_default_quadrature_order 1 -pc_type lu -dm_adaptor pragmatic
339071b71afSMatthew G. Knepley 
340071b71afSMatthew G. Knepley   test:
341edaa2528SJoe Wallwork     suffix: 2d_met3
342cbddfcd8SJoe Wallwork     args: -dm_plex_box_faces 9,9 -met 3 -dm_adaptor pragmatic -init_dm_view -adapt_dm_view \
343071b71afSMatthew 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 \
3449e3182e6SJoe Wallwork             -dm_plex_metric_target_complexity 10000.0 -dm_adaptor pragmatic
345c4762a1bSJed Brown 
346c4762a1bSJed Brown TEST*/
347