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*2a8381b2SBarry Smith static PetscErrorCode sensor(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], PetscCtx ctx) 23d71ae5a4SJacob Faibussowitsch { 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; 313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 32071b71afSMatthew G. Knepley } 33071b71afSMatthew G. Knepley 34d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 35d71ae5a4SJacob Faibussowitsch { 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 43d0609cedSBarry 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)); 49d0609cedSBarry Smith PetscOptionsEnd(); 503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51c4762a1bSJed Brown } 52c4762a1bSJed Brown 53d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 54d71ae5a4SJacob Faibussowitsch { 55c4762a1bSJed Brown PetscFunctionBegin; 569566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 579566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 589566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 599566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "DMinit")); 609566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-init_dm_view")); 613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62c4762a1bSJed Brown } 63071b71afSMatthew G. Knepley 64d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputeMetricSensor(DM dm, AppCtx *user, Vec *metric) 65d71ae5a4SJacob Faibussowitsch { 668434afd1SBarry Smith PetscSimplePointFn *funcs[1] = {sensor}; 67edaa2528SJoe Wallwork DM dmSensor, dmGrad, dmHess, dmDet; 68071b71afSMatthew G. Knepley PetscFE fe; 69edaa2528SJoe Wallwork Vec f, g, H, determinant; 70071b71afSMatthew G. Knepley PetscBool simplex; 71071b71afSMatthew G. Knepley PetscInt dim; 72c4762a1bSJed Brown 73071b71afSMatthew G. Knepley PetscFunctionBegin; 749566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 759566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 76071b71afSMatthew G. Knepley 779566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmSensor)); 789566063dSJacob Faibussowitsch PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, simplex, 1, -1, &fe)); 799566063dSJacob Faibussowitsch PetscCall(DMSetField(dmSensor, 0, NULL, (PetscObject)fe)); 809566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 819566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmSensor)); 829566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dmSensor, &f)); 839566063dSJacob Faibussowitsch PetscCall(DMProjectFunctionLocal(dmSensor, 0., funcs, NULL, INSERT_VALUES, f)); 849566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(f, NULL, "-sensor_view")); 85071b71afSMatthew G. Knepley 86071b71afSMatthew G. Knepley // Recover the gradient of the sensor function 879566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmGrad)); 889566063dSJacob Faibussowitsch PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, dim, simplex, 1, -1, &fe)); 899566063dSJacob Faibussowitsch PetscCall(DMSetField(dmGrad, 0, NULL, (PetscObject)fe)); 909566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 919566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmGrad)); 929566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(dmGrad, &g)); 939566063dSJacob Faibussowitsch PetscCall(DMPlexComputeGradientClementInterpolant(dmSensor, f, g)); 949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&f)); 959566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(g, NULL, "-gradient_view")); 96071b71afSMatthew G. Knepley 97071b71afSMatthew G. Knepley // Recover the Hessian of the sensor function 989566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmHess)); 999566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreate(dmHess, 0, &H)); 1009566063dSJacob Faibussowitsch PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, g, H)); 1019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&g)); 1029566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(H, NULL, "-hessian_view")); 103071b71afSMatthew G. Knepley 104071b71afSMatthew G. Knepley // Obtain a metric by Lp normalization 105edaa2528SJoe Wallwork PetscCall(DMPlexMetricCreate(dm, 0, metric)); 106edaa2528SJoe Wallwork PetscCall(DMPlexMetricDeterminantCreate(dm, 0, &determinant, &dmDet)); 107edaa2528SJoe Wallwork PetscCall(DMPlexMetricNormalize(dmHess, H, PETSC_TRUE, PETSC_TRUE, *metric, determinant)); 108edaa2528SJoe Wallwork PetscCall(VecDestroy(&determinant)); 109edaa2528SJoe Wallwork PetscCall(DMDestroy(&dmDet)); 1109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&H)); 1119566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmHess)); 1129566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmGrad)); 1139566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmSensor)); 1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown 117d71ae5a4SJacob Faibussowitsch static PetscErrorCode ComputeMetric(DM dm, AppCtx *user, Vec *metric) 118d71ae5a4SJacob Faibussowitsch { 119071b71afSMatthew G. Knepley PetscReal lambda = 1 / (user->hmax * user->hmax); 120c4762a1bSJed Brown 121c4762a1bSJed Brown PetscFunctionBeginUser; 122071b71afSMatthew G. Knepley if (user->metOpt == 0) { 1233e32e87aSJoe Wallwork /* Specify a uniform, isotropic metric */ 1249566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreateUniform(dm, 0, lambda, metric)); 125071b71afSMatthew G. Knepley } else if (user->metOpt == 3) { 1269566063dSJacob Faibussowitsch PetscCall(ComputeMetricSensor(dm, user, metric)); 127071b71afSMatthew G. Knepley } else { 1283e32e87aSJoe Wallwork DM cdm; 1293e32e87aSJoe Wallwork Vec coordinates; 1303e32e87aSJoe Wallwork const PetscScalar *coords; 1313e32e87aSJoe Wallwork PetscScalar *met; 1323e32e87aSJoe Wallwork PetscReal h; 133fdc00aefSJoe Wallwork PetscInt dim, i, j, vStart, vEnd, v; 1343e32e87aSJoe Wallwork 1359566063dSJacob Faibussowitsch PetscCall(DMPlexMetricCreate(dm, 0, metric)); 1369566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1379566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 1389566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 1399566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coordinates, &coords)); 1409566063dSJacob Faibussowitsch PetscCall(VecGetArray(*metric, &met)); 1419566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 142071b71afSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 143071b71afSMatthew G. Knepley PetscScalar *vcoords; 144c4762a1bSJed Brown PetscScalar *pmet; 145c4762a1bSJed Brown 1469566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRead(cdm, v, coords, &vcoords)); 147c4762a1bSJed Brown switch (user->metOpt) { 148d71ae5a4SJacob Faibussowitsch case 1: 149d71ae5a4SJacob Faibussowitsch h = user->hmax - (user->hmax - user->hmin) * PetscRealPart(vcoords[0]); 150d71ae5a4SJacob Faibussowitsch break; 151d71ae5a4SJacob Faibussowitsch case 2: 152d71ae5a4SJacob Faibussowitsch h = user->hmax * PetscAbsReal(((PetscReal)1.0) - PetscExpReal(-PetscAbsScalar(vcoords[0] - (PetscReal)0.5))) + user->hmin; 153d71ae5a4SJacob Faibussowitsch break; 154d71ae5a4SJacob Faibussowitsch default: 155d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "metOpt = 0, 1, 2 or 3, cannot be %d", user->metOpt); 156c4762a1bSJed Brown } 1579566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(dm, v, met, &pmet)); 158fdc00aefSJoe Wallwork for (i = 0; i < dim; ++i) { 159fdc00aefSJoe Wallwork for (j = 0; j < dim; ++j) { 160fdc00aefSJoe Wallwork if (i == j) { 161fdc00aefSJoe Wallwork if (i == 0) pmet[i * dim + j] = 1 / (h * h); 162fdc00aefSJoe Wallwork else pmet[i * dim + j] = lambda; 163fdc00aefSJoe Wallwork } else pmet[i * dim + j] = 0.0; 164fdc00aefSJoe Wallwork } 165fdc00aefSJoe Wallwork } 166c4762a1bSJed Brown } 1679566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(*metric, &met)); 1689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coordinates, &coords)); 1693e32e87aSJoe Wallwork } 1703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173*2a8381b2SBarry Smith static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, PetscCtx ctx) 174d71ae5a4SJacob Faibussowitsch { 175c4762a1bSJed Brown u[0] = x[0] + x[1]; 176c4762a1bSJed Brown return 0; 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestL2Projection(DM dm, DM dma, AppCtx *user) 180d71ae5a4SJacob Faibussowitsch { 181071b71afSMatthew G. Knepley PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *) = {linear}; 182071b71afSMatthew G. Knepley DM dmProj, dmaProj; 183c4762a1bSJed Brown PetscFE fe; 184071b71afSMatthew G. Knepley KSP ksp; 185071b71afSMatthew G. Knepley Mat Interp, mass, mass2; 186071b71afSMatthew G. Knepley Vec u, ua, scaling, rhs, uproj; 187c4762a1bSJed Brown PetscReal error; 188071b71afSMatthew G. Knepley PetscBool simplex; 189c4762a1bSJed Brown PetscInt dim; 190c4762a1bSJed Brown 191c4762a1bSJed Brown PetscFunctionBeginUser; 1929566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 1939566063dSJacob Faibussowitsch PetscCall(DMPlexIsSimplex(dm, &simplex)); 194c4762a1bSJed Brown 1959566063dSJacob Faibussowitsch PetscCall(DMClone(dm, &dmProj)); 1969566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe)); 1979566063dSJacob Faibussowitsch PetscCall(DMSetField(dmProj, 0, NULL, (PetscObject)fe)); 1989566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 1999566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmProj)); 200071b71afSMatthew G. Knepley 2019566063dSJacob Faibussowitsch PetscCall(DMClone(dma, &dmaProj)); 2029566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe)); 2039566063dSJacob Faibussowitsch PetscCall(DMSetField(dmaProj, 0, NULL, (PetscObject)fe)); 2049566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 2059566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dmaProj)); 206071b71afSMatthew G. Knepley 2079566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmProj, &u)); 2089566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmaProj, &ua)); 2099566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmaProj, &rhs)); 2109566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dmaProj, &uproj)); 211071b71afSMatthew G. Knepley 212071b71afSMatthew G. Knepley // Interpolate onto original mesh using dual basis 2139566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dmProj, 0.0, funcs, NULL, INSERT_VALUES, u)); 2149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, "Original")); 2159566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(u, NULL, "-orig_vec_view")); 2169566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dmProj, 0.0, funcs, NULL, u, &error)); 2179566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Original L2 Error: %g\n", (double)error)); 218071b71afSMatthew G. Knepley // Interpolate onto NEW mesh using dual basis 2199566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dmaProj, 0.0, funcs, NULL, INSERT_VALUES, ua)); 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)ua, "Adapted")); 2219566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ua, NULL, "-adapt_vec_view")); 2229566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, ua, &error)); 2239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Adapted L2 Error: %g\n", (double)error)); 224071b71afSMatthew G. Knepley // Interpolate between meshes using interpolation matrix 2259566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(dmProj, dmaProj, &Interp, &scaling)); 2269566063dSJacob Faibussowitsch PetscCall(MatInterpolate(Interp, u, ua)); 2279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Interp)); 2289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&scaling)); 2299566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)ua, "Interpolation")); 2309566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ua, NULL, "-interp_vec_view")); 2319566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, ua, &error)); 2329566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Interpolated L2 Error: %g\n", (double)error)); 233071b71afSMatthew G. Knepley // L2 projection 2349566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dmaProj, dmaProj, &mass)); 2359566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(mass, NULL, "-mass_mat_view")); 2369566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp)); 2379566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, mass, mass)); 2389566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 239da81f932SPierre Jolivet // Compute rhs as M f, could also directly project the analytic function but we might not have it 2409566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dmProj, dmaProj, &mass2)); 2419566063dSJacob Faibussowitsch PetscCall(MatMult(mass2, u, rhs)); 2429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mass2)); 2439566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, rhs, uproj)); 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)uproj, "L_2 Projection")); 2459566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(uproj, NULL, "-proj_vec_view")); 2469566063dSJacob Faibussowitsch PetscCall(DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, uproj, &error)); 2479566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double)error)); 2489566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 2499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mass)); 2509566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmProj, &u)); 2519566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmaProj, &ua)); 2529566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmaProj, &rhs)); 2539566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dmaProj, &uproj)); 2549566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmProj)); 2559566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmaProj)); 2563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 257c4762a1bSJed Brown } 258c4762a1bSJed Brown 259d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[]) 260d71ae5a4SJacob Faibussowitsch { 261071b71afSMatthew G. Knepley DM dm; 262c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 263c4762a1bSJed Brown MPI_Comm comm; 264c4762a1bSJed Brown DM dma, odm; 265c4762a1bSJed Brown Vec metric; 266071b71afSMatthew G. Knepley PetscInt r; 267c4762a1bSJed Brown 268327415f7SBarry Smith PetscFunctionBeginUser; 2699566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 270c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 2719566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user)); 2729566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &dm)); 273c4762a1bSJed Brown 274071b71afSMatthew G. Knepley odm = dm; 2759566063dSJacob Faibussowitsch PetscCall(DMPlexDistributeOverlap(odm, 1, NULL, &dm)); 2769371c9d4SSatish Balay if (!dm) { 2779371c9d4SSatish Balay dm = odm; 2789371c9d4SSatish Balay } else PetscCall(DMDestroy(&odm)); 279071b71afSMatthew G. Knepley 280071b71afSMatthew G. Knepley for (r = 0; r < user.Nr; ++r) { 281071b71afSMatthew G. Knepley DMLabel label; 282071b71afSMatthew G. Knepley 2839566063dSJacob Faibussowitsch PetscCall(ComputeMetric(dm, &user, &metric)); 2849566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "marker", &label)); 2859566063dSJacob Faibussowitsch PetscCall(DMAdaptMetric(dm, metric, label, NULL, &dma)); 2869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&metric)); 2879566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dma, "DMadapt")); 2889566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dma, "adapt_")); 2899566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dma, NULL, "-dm_view")); 2909566063dSJacob Faibussowitsch if (user.doL2) PetscCall(TestL2Projection(dm, dma, &user)); 2919566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 292071b71afSMatthew G. Knepley dm = dma; 293071b71afSMatthew G. Knepley } 2949566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, "final_")); 2959566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 2969566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 2979566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 298b122ec5aSJacob Faibussowitsch return 0; 299c4762a1bSJed Brown } 300c4762a1bSJed Brown 301c4762a1bSJed Brown /*TEST 302c4762a1bSJed Brown 3038e3a2eefSMatthew G. Knepley build: 3048e3a2eefSMatthew G. Knepley requires: pragmatic 3058e3a2eefSMatthew G. Knepley 306071b71afSMatthew G. Knepley testset: 3079e3182e6SJoe Wallwork args: -dm_plex_box_faces 4,4,4 -dm_adaptor pragmatic -met 2 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic 308071b71afSMatthew G. Knepley 309c4762a1bSJed Brown test: 310edaa2528SJoe Wallwork suffix: 2d 311071b71afSMatthew G. Knepley args: -dm_plex_separate_marker 0 312c4762a1bSJed Brown test: 313edaa2528SJoe Wallwork suffix: 2d_sep 314071b71afSMatthew G. Knepley args: -dm_plex_separate_marker 1 315c4762a1bSJed Brown test: 316edaa2528SJoe Wallwork suffix: 3d 317071b71afSMatthew G. Knepley args: -dm_plex_dim 3 318071b71afSMatthew G. Knepley 319071b71afSMatthew G. Knepley # Pragmatic hangs for simple partitioner 320071b71afSMatthew G. Knepley testset: 321071b71afSMatthew G. Knepley requires: parmetis 322edaa2528SJoe Wallwork args: -dm_plex_box_faces 2,2 -petscpartitioner_type parmetis -met 2 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic 323071b71afSMatthew G. Knepley 324c4762a1bSJed Brown test: 325edaa2528SJoe Wallwork suffix: 2d_parmetis_np2 326c4762a1bSJed Brown nsize: 2 327c4762a1bSJed Brown test: 328edaa2528SJoe Wallwork suffix: 2d_parmetis_np4 329c4762a1bSJed Brown nsize: 4 330071b71afSMatthew G. Knepley 331c4762a1bSJed Brown test: 3329880b877SBarry Smith requires: parmetis 333edaa2528SJoe Wallwork suffix: 3d_parmetis_met0 334c4762a1bSJed Brown nsize: 2 335e600fa54SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 9,9,9 -dm_adaptor pragmatic -petscpartitioner_type parmetis \ 3369e3182e6SJoe Wallwork -met 0 -hmin 0.01 -hmax 0.03 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic 337c4762a1bSJed Brown test: 3389880b877SBarry Smith requires: parmetis 339edaa2528SJoe Wallwork suffix: 3d_parmetis_met2 3403e32e87aSJoe Wallwork nsize: 2 341e600fa54SMatthew G. Knepley args: -dm_plex_box_faces 19,19 -dm_adaptor pragmatic -petscpartitioner_type parmetis \ 3429e3182e6SJoe Wallwork -met 2 -hmax 0.5 -hmin 0.001 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic 343c4762a1bSJed Brown test: 344edaa2528SJoe Wallwork suffix: proj2 345cbddfcd8SJoe Wallwork args: -dm_plex_box_faces 2,2 -dm_plex_hash_location -dm_adaptor pragmatic -init_dm_view -adapt_dm_view -do_L2 \ 3469e3182e6SJoe Wallwork -petscspace_degree 1 -petscfe_default_quadrature_order 1 -pc_type lu -dm_adaptor pragmatic 347c4762a1bSJed Brown test: 348edaa2528SJoe Wallwork suffix: proj4 349cbddfcd8SJoe Wallwork args: -dm_plex_box_faces 4,4 -dm_plex_hash_location -dm_adaptor pragmatic -init_dm_view -adapt_dm_view -do_L2 \ 3509e3182e6SJoe Wallwork -petscspace_degree 1 -petscfe_default_quadrature_order 1 -pc_type lu -dm_adaptor pragmatic 351071b71afSMatthew G. Knepley 352071b71afSMatthew G. Knepley test: 353edaa2528SJoe Wallwork suffix: 2d_met3 354cbddfcd8SJoe Wallwork args: -dm_plex_box_faces 9,9 -met 3 -dm_adaptor pragmatic -init_dm_view -adapt_dm_view \ 355071b71afSMatthew 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 \ 3569e3182e6SJoe Wallwork -dm_plex_metric_target_complexity 10000.0 -dm_adaptor pragmatic 357c4762a1bSJed Brown 358c4762a1bSJed Brown TEST*/ 359