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 PetscFunctionBegin; 37071b71afSMatthew G. Knepley options->Nr = 1; 38c4762a1bSJed Brown options->metOpt = 1; 39c4762a1bSJed Brown options->hmin = 0.05; 40c4762a1bSJed Brown options->hmax = 0.5; 41c4762a1bSJed Brown options->doL2 = PETSC_FALSE; 42c4762a1bSJed Brown 43*d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Adaptation Options", "DMPLEX"); 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-Nr", "Numberof refinement passes", "ex19.c", options->Nr, &options->Nr, NULL, 1)); 459566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-met", "Different choices of metric", "ex19.c", options->metOpt, &options->metOpt, NULL,0)); 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-hmax", "Max size prescribed by the metric", "ex19.c", options->hmax, &options->hmax, NULL)); 479566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-hmin", "Min size prescribed by the metric", "ex19.c", options->hmin, &options->hmin, NULL)); 489566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-do_L2", "Test L2 projection", "ex19.c", options->doL2, &options->doL2, NULL)); 49*d0609cedSBarry Smith PetscOptionsEnd(); 50c4762a1bSJed Brown 51c4762a1bSJed Brown PetscFunctionReturn(0); 52c4762a1bSJed Brown } 53c4762a1bSJed Brown 54071b71afSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 55c4762a1bSJed Brown { 56c4762a1bSJed Brown PetscFunctionBegin; 579566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 589566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 599566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *dm, "DMinit")); 619566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-init_dm_view")); 62071b71afSMatthew G. Knepley PetscFunctionReturn(0); 63c4762a1bSJed Brown } 64071b71afSMatthew G. Knepley 65071b71afSMatthew G. Knepley static PetscErrorCode ComputeMetricSensor(DM dm, AppCtx *user, Vec *metric) 66c4762a1bSJed Brown { 67071b71afSMatthew G. Knepley PetscSimplePointFunc funcs[1] = {sensor}; 68071b71afSMatthew G. Knepley DM dmSensor, dmGrad, dmHess; 69071b71afSMatthew G. Knepley PetscFE fe; 70071b71afSMatthew G. Knepley Vec f, g, H; 71071b71afSMatthew G. Knepley PetscBool simplex; 72071b71afSMatthew G. Knepley PetscInt dim; 73c4762a1bSJed Brown 74071b71afSMatthew G. Knepley PetscFunctionBegin; 759566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 769566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 77071b71afSMatthew G. Knepley 789566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmSensor)); 799566063dSJacob Faibussowitsch PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, simplex, 1, -1, &fe)); 809566063dSJacob Faibussowitsch PetscCall(DMSetField(dmSensor, 0, NULL, (PetscObject) fe)); 819566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 829566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmSensor)); 839566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dmSensor, &f)); 849566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dmSensor, 0., funcs, NULL, INSERT_VALUES, f)); 859566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(f, NULL, "-sensor_view")); 86071b71afSMatthew G. Knepley 87071b71afSMatthew G. Knepley // Recover the gradient of the sensor function 889566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmGrad)); 899566063dSJacob Faibussowitsch PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, dim, simplex, 1, -1, &fe)); 909566063dSJacob Faibussowitsch PetscCall(DMSetField(dmGrad, 0, NULL, (PetscObject) fe)); 919566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 929566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmGrad)); 939566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dmGrad, &g)); 949566063dSJacob Faibussowitsch PetscCall(DMPlexComputeGradientClementInterpolant(dmSensor, f, g)); 959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&f)); 969566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(g, NULL, "-gradient_view")); 97071b71afSMatthew G. Knepley 98071b71afSMatthew G. Knepley // Recover the Hessian of the sensor function 999566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmHess)); 1009566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreate(dmHess, 0, &H)); 1019566063dSJacob Faibussowitsch PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, g, H)); 1029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g)); 1039566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(H, NULL, "-hessian_view")); 104071b71afSMatthew G. Knepley 105071b71afSMatthew G. Knepley // Obtain a metric by Lp normalization 1069566063dSJacob Faibussowitsch PetscCall(DMPlexMetricNormalize(dmHess, H, PETSC_TRUE, PETSC_TRUE, metric)); 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 114c4762a1bSJed Brown static PetscErrorCode ComputeMetric(DM dm, AppCtx *user, Vec *metric) 115c4762a1bSJed Brown { 116071b71afSMatthew G. Knepley PetscReal lambda = 1/(user->hmax*user->hmax); 117c4762a1bSJed Brown 118c4762a1bSJed Brown PetscFunctionBeginUser; 119071b71afSMatthew G. Knepley if (user->metOpt == 0) { 1203e32e87aSJoe Wallwork /* Specify a uniform, isotropic metric */ 1219566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreateUniform(dm, 0, lambda, metric)); 122071b71afSMatthew G. Knepley } else if (user->metOpt == 3) { 1239566063dSJacob Faibussowitsch PetscCall(ComputeMetricSensor(dm, user, metric)); 124071b71afSMatthew G. Knepley } else { 1253e32e87aSJoe Wallwork DM cdm; 1263e32e87aSJoe Wallwork Vec coordinates; 1273e32e87aSJoe Wallwork const PetscScalar *coords; 1283e32e87aSJoe Wallwork PetscScalar *met; 1293e32e87aSJoe Wallwork PetscReal h; 130fdc00aefSJoe Wallwork PetscInt dim, i, j, vStart, vEnd, v; 1313e32e87aSJoe Wallwork 1329566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreate(dm, 0, metric)); 1339566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1349566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 1359566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordinates, &coords)); 1379566063dSJacob Faibussowitsch PetscCall(VecGetArray(*metric, &met)); 1389566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 139071b71afSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 140071b71afSMatthew G. Knepley PetscScalar *vcoords; 141c4762a1bSJed Brown PetscScalar *pmet; 142c4762a1bSJed Brown 1439566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(cdm, v, coords, &vcoords)); 144c4762a1bSJed Brown switch (user->metOpt) { 145c4762a1bSJed Brown case 1: 146071b71afSMatthew G. Knepley h = user->hmax - (user->hmax-user->hmin)*PetscRealPart(vcoords[0]); 147c4762a1bSJed Brown break; 148c4762a1bSJed Brown case 2: 149071b71afSMatthew G. Knepley h = user->hmax*PetscAbsReal(((PetscReal) 1.0)-PetscExpReal(-PetscAbsScalar(vcoords[0]-(PetscReal)0.5))) + user->hmin; 150c4762a1bSJed Brown break; 151c4762a1bSJed Brown default: 15298921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "metOpt = 0, 1, 2 or 3, cannot be %d", user->metOpt); 153c4762a1bSJed Brown } 1549566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dm, v, met, &pmet)); 155fdc00aefSJoe Wallwork for (i = 0; i < dim; ++i) { 156fdc00aefSJoe Wallwork for (j = 0; j < dim; ++j) { 157fdc00aefSJoe Wallwork if (i == j) { 158fdc00aefSJoe Wallwork if (i == 0) pmet[i*dim+j] = 1/(h*h); 159fdc00aefSJoe Wallwork else pmet[i*dim+j] = lambda; 160fdc00aefSJoe Wallwork } else pmet[i*dim+j] = 0.0; 161fdc00aefSJoe Wallwork } 162fdc00aefSJoe Wallwork } 163c4762a1bSJed Brown } 1649566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*metric, &met)); 1659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordinates, &coords)); 1663e32e87aSJoe Wallwork } 167c4762a1bSJed Brown PetscFunctionReturn(0); 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170c4762a1bSJed Brown static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 171c4762a1bSJed Brown { 172c4762a1bSJed Brown u[0] = x[0] + x[1]; 173c4762a1bSJed Brown return 0; 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown static PetscErrorCode TestL2Projection(DM dm, DM dma, AppCtx *user) 177c4762a1bSJed Brown { 178071b71afSMatthew G. Knepley PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = {linear}; 179071b71afSMatthew G. Knepley DM dmProj, dmaProj; 180c4762a1bSJed Brown PetscFE fe; 181071b71afSMatthew G. Knepley KSP ksp; 182071b71afSMatthew G. Knepley Mat Interp, mass, mass2; 183071b71afSMatthew G. Knepley Vec u, ua, scaling, rhs, uproj; 184c4762a1bSJed Brown PetscReal error; 185071b71afSMatthew G. Knepley PetscBool simplex; 186c4762a1bSJed Brown PetscInt dim; 187c4762a1bSJed Brown 188c4762a1bSJed Brown PetscFunctionBeginUser; 1899566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1909566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 191c4762a1bSJed Brown 1929566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmProj)); 1939566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe)); 1949566063dSJacob Faibussowitsch PetscCall(DMSetField(dmProj, 0, NULL, (PetscObject) fe)); 1959566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 1969566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmProj)); 197071b71afSMatthew G. Knepley 1989566063dSJacob Faibussowitsch PetscCall(DMClone(dma, &dmaProj)); 1999566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe)); 2009566063dSJacob Faibussowitsch PetscCall(DMSetField(dmaProj, 0, NULL, (PetscObject) fe)); 2019566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 2029566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmaProj)); 203071b71afSMatthew G. Knepley 2049566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmProj, &u)); 2059566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmaProj, &ua)); 2069566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmaProj, &rhs)); 2079566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmaProj, &uproj)); 208071b71afSMatthew G. Knepley 209071b71afSMatthew G. Knepley // Interpolate onto original mesh using dual basis 2109566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dmProj, 0.0, funcs, NULL, INSERT_VALUES, u)); 2119566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) u, "Original")); 2129566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-orig_vec_view")); 2139566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dmProj, 0.0, funcs, NULL, u, &error)); 2149566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Original L2 Error: %g\n", (double) error)); 215071b71afSMatthew G. Knepley // Interpolate onto NEW mesh using dual basis 2169566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dmaProj, 0.0, funcs, NULL, INSERT_VALUES, ua)); 2179566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) ua, "Adapted")); 2189566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ua, NULL, "-adapt_vec_view")); 2199566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, ua, &error)); 2209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Adapted L2 Error: %g\n", (double) error)); 221071b71afSMatthew G. Knepley // Interpolate between meshes using interpolation matrix 2229566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dmProj, dmaProj, &Interp, &scaling)); 2239566063dSJacob Faibussowitsch PetscCall(MatInterpolate(Interp, u, ua)); 2249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Interp)); 2259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scaling)); 2269566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) ua, "Interpolation")); 2279566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ua, NULL, "-interp_vec_view")); 2289566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, ua, &error)); 2299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Interpolated L2 Error: %g\n", (double) error)); 230071b71afSMatthew G. Knepley // L2 projection 2319566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dmaProj, dmaProj, &mass)); 2329566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(mass, NULL, "-mass_mat_view")); 2339566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp)); 2349566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, mass, mass)); 2359566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 236071b71afSMatthew G. Knepley // Compute rhs as M f, could also direclty project the analytic function but we might not have it 2379566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dmProj, dmaProj, &mass2)); 2389566063dSJacob Faibussowitsch PetscCall(MatMult(mass2, u, rhs)); 2399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mass2)); 2409566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, rhs, uproj)); 2419566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) uproj, "L_2 Projection")); 2429566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uproj, NULL, "-proj_vec_view")); 2439566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, uproj, &error)); 2449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double) error)); 2459566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 2469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mass)); 2479566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmProj, &u)); 2489566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmaProj, &ua)); 2499566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmaProj, &rhs)); 2509566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmaProj, &uproj)); 2519566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmProj)); 2529566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmaProj)); 253c4762a1bSJed Brown PetscFunctionReturn(0); 254c4762a1bSJed Brown } 255c4762a1bSJed Brown 256c4762a1bSJed Brown int main (int argc, char * argv[]) { 257071b71afSMatthew G. Knepley DM dm; 258c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 259c4762a1bSJed Brown MPI_Comm comm; 260c4762a1bSJed Brown DM dma, odm; 261c4762a1bSJed Brown Vec metric; 262071b71afSMatthew G. Knepley PetscInt r; 263c4762a1bSJed Brown 2649566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 265c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 2669566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 2679566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &dm)); 268c4762a1bSJed Brown 269071b71afSMatthew G. Knepley odm = dm; 2709566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeOverlap(odm, 1, NULL, &dm)); 271071b71afSMatthew G. Knepley if (!dm) {dm = odm;} 2729566063dSJacob Faibussowitsch else PetscCall(DMDestroy(&odm)); 273071b71afSMatthew G. Knepley 274071b71afSMatthew G. Knepley for (r = 0; r < user.Nr; ++r) { 275071b71afSMatthew G. Knepley DMLabel label; 276071b71afSMatthew G. Knepley 2779566063dSJacob Faibussowitsch PetscCall(ComputeMetric(dm, &user, &metric)); 2789566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 2799566063dSJacob Faibussowitsch PetscCall(DMAdaptMetric(dm, metric, label, NULL, &dma)); 2809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metric)); 2819566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) dma, "DMadapt")); 2829566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) dma, "adapt_")); 2839566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dma, NULL, "-dm_view")); 2849566063dSJacob Faibussowitsch if (user.doL2) PetscCall(TestL2Projection(dm, dma, &user)); 2859566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 286071b71afSMatthew G. Knepley dm = dma; 287071b71afSMatthew G. Knepley } 2889566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) dm, "final_")); 2899566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 2909566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2919566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 292b122ec5aSJacob Faibussowitsch return 0; 293c4762a1bSJed Brown } 294c4762a1bSJed Brown 295c4762a1bSJed Brown /*TEST 296c4762a1bSJed Brown 2978e3a2eefSMatthew G. Knepley build: 2988e3a2eefSMatthew G. Knepley requires: pragmatic 2998e3a2eefSMatthew G. Knepley 300071b71afSMatthew G. Knepley testset: 3019e3182e6SJoe Wallwork args: -dm_plex_box_faces 4,4,4 -dm_adaptor pragmatic -met 2 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic 302071b71afSMatthew G. Knepley 303c4762a1bSJed Brown test: 304c4762a1bSJed Brown suffix: 0 305071b71afSMatthew G. Knepley args: -dm_plex_separate_marker 0 306c4762a1bSJed Brown test: 307c4762a1bSJed Brown suffix: 1 308071b71afSMatthew G. Knepley args: -dm_plex_separate_marker 1 309c4762a1bSJed Brown test: 310c4762a1bSJed Brown suffix: 2 311071b71afSMatthew G. Knepley args: -dm_plex_dim 3 312c4762a1bSJed Brown test: 313c4762a1bSJed Brown suffix: 3 314071b71afSMatthew G. Knepley args: -dm_plex_dim 3 315071b71afSMatthew G. Knepley 316071b71afSMatthew G. Knepley # Pragmatic hangs for simple partitioner 317071b71afSMatthew G. Knepley testset: 318071b71afSMatthew G. Knepley requires: parmetis 319e600fa54SMatthew 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 320071b71afSMatthew G. Knepley 321c4762a1bSJed Brown test: 322c4762a1bSJed Brown suffix: 4 323c4762a1bSJed Brown nsize: 2 324c4762a1bSJed Brown test: 325c4762a1bSJed Brown suffix: 5 326c4762a1bSJed Brown nsize: 4 327071b71afSMatthew G. Knepley 328c4762a1bSJed Brown test: 3299880b877SBarry Smith requires: parmetis 330c4762a1bSJed Brown suffix: 6 331c4762a1bSJed Brown nsize: 2 332e600fa54SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 9,9,9 -dm_adaptor pragmatic -petscpartitioner_type parmetis \ 3339e3182e6SJoe Wallwork -met 0 -hmin 0.01 -hmax 0.03 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic 334c4762a1bSJed Brown test: 3359880b877SBarry Smith requires: parmetis 336c4762a1bSJed Brown suffix: 7 3373e32e87aSJoe Wallwork nsize: 2 338e600fa54SMatthew G. Knepley args: -dm_plex_box_faces 19,19 -dm_adaptor pragmatic -petscpartitioner_type parmetis \ 3399e3182e6SJoe Wallwork -met 2 -hmax 0.5 -hmin 0.001 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic 340c4762a1bSJed Brown test: 341c4762a1bSJed Brown suffix: proj_0 342cbddfcd8SJoe Wallwork args: -dm_plex_box_faces 2,2 -dm_plex_hash_location -dm_adaptor pragmatic -init_dm_view -adapt_dm_view -do_L2 \ 3439e3182e6SJoe Wallwork -petscspace_degree 1 -petscfe_default_quadrature_order 1 -pc_type lu -dm_adaptor pragmatic 344c4762a1bSJed Brown test: 345c4762a1bSJed Brown suffix: proj_1 346cbddfcd8SJoe Wallwork args: -dm_plex_box_faces 4,4 -dm_plex_hash_location -dm_adaptor pragmatic -init_dm_view -adapt_dm_view -do_L2 \ 3479e3182e6SJoe Wallwork -petscspace_degree 1 -petscfe_default_quadrature_order 1 -pc_type lu -dm_adaptor pragmatic 348071b71afSMatthew G. Knepley 349071b71afSMatthew G. Knepley test: 350071b71afSMatthew G. Knepley suffix: sensor 351cbddfcd8SJoe Wallwork args: -dm_plex_box_faces 9,9 -met 3 -dm_adaptor pragmatic -init_dm_view -adapt_dm_view \ 352071b71afSMatthew 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 \ 3539e3182e6SJoe Wallwork -dm_plex_metric_target_complexity 10000.0 -dm_adaptor pragmatic 354c4762a1bSJed Brown 355c4762a1bSJed Brown TEST*/ 356