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 { 8*071b71afSMatthew 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 14*071b71afSMatthew G. Knepley /* 15*071b71afSMatthew G. Knepley Classic hyperbolic sensor function for testing multi-scale anisotropic mesh adaptation: 16*071b71afSMatthew G. Knepley 17*071b71afSMatthew G. Knepley f:[-1, 1]x[-1, 1] \to R, 18*071b71afSMatthew G. Knepley f(x, y) = sin(50xy)/100 if |xy| > 2\pi/50 else sin(50xy) 19*071b71afSMatthew G. Knepley 20*071b71afSMatthew G. Knepley (mapped to have domain [0,1] x [0,1] in this case). 21*071b71afSMatthew G. Knepley */ 22*071b71afSMatthew G. Knepley static PetscErrorCode sensor(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar u[], void *ctx) 23*071b71afSMatthew G. Knepley { 24*071b71afSMatthew G. Knepley const PetscReal xref = 2.*x[0] - 1.; 25*071b71afSMatthew G. Knepley const PetscReal yref = 2.*x[1] - 1.; 26*071b71afSMatthew G. Knepley const PetscReal xy = xref*yref; 27*071b71afSMatthew G. Knepley 28*071b71afSMatthew G. Knepley PetscFunctionBeginUser; 29*071b71afSMatthew G. Knepley u[0] = PetscSinReal(50.*xy); 30*071b71afSMatthew G. Knepley if (PetscAbsReal(xy) > 2.*PETSC_PI/50.) u[0] *= 0.01; 31*071b71afSMatthew G. Knepley PetscFunctionReturn(0); 32*071b71afSMatthew G. Knepley } 33*071b71afSMatthew G. Knepley 34c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 35c4762a1bSJed Brown { 36c4762a1bSJed Brown PetscErrorCode ierr; 37c4762a1bSJed Brown 38c4762a1bSJed Brown PetscFunctionBegin; 39*071b71afSMatthew 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 45c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Meshing Adaptation Options", "DMPLEX");CHKERRQ(ierr); 46*071b71afSMatthew G. Knepley ierr = PetscOptionsBoundedInt("-Nr", "Numberof refinement passes", "ex19.c", options->Nr, &options->Nr, NULL, 1);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-met", "Different choices of metric", "ex19.c", options->metOpt, &options->metOpt, NULL,0);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = PetscOptionsReal("-hmax", "Max size prescribed by the metric", "ex19.c", options->hmax, &options->hmax, NULL);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = PetscOptionsReal("-hmin", "Min size prescribed by the metric", "ex19.c", options->hmin, &options->hmin, NULL);CHKERRQ(ierr); 50c4762a1bSJed Brown ierr = PetscOptionsBool("-do_L2", "Test L2 projection", "ex19.c", options->doL2, &options->doL2, NULL);CHKERRQ(ierr); 511e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 52c4762a1bSJed Brown 53c4762a1bSJed Brown PetscFunctionReturn(0); 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56*071b71afSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 57c4762a1bSJed Brown { 58c4762a1bSJed Brown PetscErrorCode ierr; 59c4762a1bSJed Brown 60c4762a1bSJed Brown PetscFunctionBegin; 61*071b71afSMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 62*071b71afSMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 63*071b71afSMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 64*071b71afSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *dm, "DMinit");CHKERRQ(ierr); 65*071b71afSMatthew G. Knepley ierr = DMViewFromOptions(*dm, NULL, "-init_dm_view");CHKERRQ(ierr); 66*071b71afSMatthew G. Knepley PetscFunctionReturn(0); 67c4762a1bSJed Brown } 68*071b71afSMatthew G. Knepley 69*071b71afSMatthew G. Knepley static PetscErrorCode ComputeMetricSensor(DM dm, AppCtx *user, Vec *metric) 70c4762a1bSJed Brown { 71*071b71afSMatthew G. Knepley PetscSimplePointFunc funcs[1] = {sensor}; 72*071b71afSMatthew G. Knepley DM dmSensor, dmGrad, dmHess; 73*071b71afSMatthew G. Knepley PetscFE fe; 74*071b71afSMatthew G. Knepley Vec f, g, H; 75*071b71afSMatthew G. Knepley PetscBool simplex; 76*071b71afSMatthew G. Knepley PetscInt dim; 77*071b71afSMatthew G. Knepley PetscErrorCode ierr; 78c4762a1bSJed Brown 79*071b71afSMatthew G. Knepley PetscFunctionBegin; 80*071b71afSMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 81*071b71afSMatthew G. Knepley ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 82*071b71afSMatthew G. Knepley 83*071b71afSMatthew G. Knepley ierr = DMClone(dm, &dmSensor);CHKERRQ(ierr); 84*071b71afSMatthew G. Knepley ierr = PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, simplex, 1, -1, &fe);CHKERRQ(ierr); 85*071b71afSMatthew G. Knepley ierr = DMSetField(dmSensor, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 86*071b71afSMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 87*071b71afSMatthew G. Knepley ierr = DMCreateDS(dmSensor);CHKERRQ(ierr); 88*071b71afSMatthew G. Knepley ierr = DMCreateLocalVector(dmSensor, &f);CHKERRQ(ierr); 89*071b71afSMatthew G. Knepley ierr = DMProjectFunctionLocal(dmSensor, 0., funcs, NULL, INSERT_VALUES, f);CHKERRQ(ierr); 90*071b71afSMatthew G. Knepley ierr = VecViewFromOptions(f, NULL, "-sensor_view");CHKERRQ(ierr); 91*071b71afSMatthew G. Knepley 92*071b71afSMatthew G. Knepley // Recover the gradient of the sensor function 93*071b71afSMatthew G. Knepley ierr = DMClone(dm, &dmGrad);CHKERRQ(ierr); 94*071b71afSMatthew G. Knepley ierr = PetscFECreateLagrange(PETSC_COMM_SELF, dim, dim, simplex, 1, -1, &fe);CHKERRQ(ierr); 95*071b71afSMatthew G. Knepley ierr = DMSetField(dmGrad, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 96*071b71afSMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 97*071b71afSMatthew G. Knepley ierr = DMCreateDS(dmGrad);CHKERRQ(ierr); 98*071b71afSMatthew G. Knepley ierr = DMCreateLocalVector(dmGrad, &g);CHKERRQ(ierr); 99*071b71afSMatthew G. Knepley ierr = DMPlexComputeGradientClementInterpolant(dmSensor, f, g);CHKERRQ(ierr); 100*071b71afSMatthew G. Knepley ierr = VecDestroy(&f);CHKERRQ(ierr); 101*071b71afSMatthew G. Knepley ierr = VecViewFromOptions(g, NULL, "-gradient_view");CHKERRQ(ierr); 102*071b71afSMatthew G. Knepley 103*071b71afSMatthew G. Knepley // Recover the Hessian of the sensor function 104*071b71afSMatthew G. Knepley ierr = DMClone(dm, &dmHess);CHKERRQ(ierr); 105*071b71afSMatthew G. Knepley ierr = DMPlexMetricCreate(dmHess, 0, &H);CHKERRQ(ierr); 106*071b71afSMatthew G. Knepley ierr = DMPlexComputeGradientClementInterpolant(dmGrad, g, H);CHKERRQ(ierr); 107*071b71afSMatthew G. Knepley ierr = VecDestroy(&g);CHKERRQ(ierr); 108*071b71afSMatthew G. Knepley ierr = VecViewFromOptions(H, NULL, "-hessian_view");CHKERRQ(ierr); 109*071b71afSMatthew G. Knepley 110*071b71afSMatthew G. Knepley // Obtain a metric by Lp normalization 111*071b71afSMatthew G. Knepley ierr = DMPlexMetricNormalize(dmHess, H, PETSC_TRUE, PETSC_TRUE, metric);CHKERRQ(ierr); 112*071b71afSMatthew G. Knepley ierr = VecDestroy(&H);CHKERRQ(ierr); 113*071b71afSMatthew G. Knepley ierr = DMDestroy(&dmHess);CHKERRQ(ierr); 114*071b71afSMatthew G. Knepley ierr = DMDestroy(&dmGrad);CHKERRQ(ierr); 115*071b71afSMatthew G. Knepley ierr = DMDestroy(&dmSensor);CHKERRQ(ierr); 116c4762a1bSJed Brown PetscFunctionReturn(0); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown static PetscErrorCode ComputeMetric(DM dm, AppCtx *user, Vec *metric) 120c4762a1bSJed Brown { 121*071b71afSMatthew G. Knepley PetscReal lambda = 1/(user->hmax*user->hmax); 122c4762a1bSJed Brown PetscErrorCode ierr; 123c4762a1bSJed Brown 124c4762a1bSJed Brown PetscFunctionBeginUser; 125*071b71afSMatthew G. Knepley if (user->metOpt == 0) { 1263e32e87aSJoe Wallwork /* Specify a uniform, isotropic metric */ 1273e32e87aSJoe Wallwork ierr = DMPlexMetricCreateUniform(dm, 0, lambda, metric);CHKERRQ(ierr); 128*071b71afSMatthew G. Knepley } else if (user->metOpt == 3) { 129*071b71afSMatthew G. Knepley ierr = ComputeMetricSensor(dm, user, metric);CHKERRQ(ierr); 130*071b71afSMatthew G. Knepley } else { 1313e32e87aSJoe Wallwork DM cdm; 1323e32e87aSJoe Wallwork Vec coordinates; 1333e32e87aSJoe Wallwork const PetscScalar *coords; 1343e32e87aSJoe Wallwork PetscScalar *met; 1353e32e87aSJoe Wallwork PetscReal h; 136*071b71afSMatthew G. Knepley PetscInt vStart, vEnd, v; 1373e32e87aSJoe Wallwork 138*071b71afSMatthew G. Knepley ierr = DMPlexMetricCreateUniform(dm, 0, lambda, metric);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = VecGetArray(*metric, &met);CHKERRQ(ierr); 143*071b71afSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 144*071b71afSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 145*071b71afSMatthew G. Knepley PetscScalar *vcoords; 146c4762a1bSJed Brown PetscScalar *pmet; 147c4762a1bSJed Brown 148*071b71afSMatthew G. Knepley ierr = DMPlexPointLocalRead(cdm, v, coords, &vcoords);CHKERRQ(ierr); 149c4762a1bSJed Brown switch (user->metOpt) { 150c4762a1bSJed Brown case 1: 151*071b71afSMatthew G. Knepley h = user->hmax - (user->hmax-user->hmin)*PetscRealPart(vcoords[0]); 152c4762a1bSJed Brown break; 153c4762a1bSJed Brown case 2: 154*071b71afSMatthew G. Knepley h = user->hmax*PetscAbsReal(((PetscReal) 1.0)-PetscExpReal(-PetscAbsScalar(vcoords[0]-(PetscReal)0.5))) + user->hmin; 155c4762a1bSJed Brown break; 156c4762a1bSJed Brown default: 157c4762a1bSJed Brown SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "metOpt = 0, 1 or 2, cannot be %d", user->metOpt); 158c4762a1bSJed Brown } 1593e32e87aSJoe Wallwork lambda = 1/(h*h); 160*071b71afSMatthew G. Knepley ierr = DMPlexPointLocalRef(dm, v, met, &pmet);CHKERRQ(ierr); 1613e32e87aSJoe Wallwork pmet[0] = lambda; 162c4762a1bSJed Brown } 163c4762a1bSJed Brown ierr = VecRestoreArray(*metric, &met);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 1653e32e87aSJoe Wallwork } 166c4762a1bSJed Brown PetscFunctionReturn(0); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 169c4762a1bSJed Brown static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 170c4762a1bSJed Brown { 171c4762a1bSJed Brown u[0] = x[0] + x[1]; 172c4762a1bSJed Brown return 0; 173c4762a1bSJed Brown } 174c4762a1bSJed Brown 175c4762a1bSJed Brown static PetscErrorCode TestL2Projection(DM dm, DM dma, AppCtx *user) 176c4762a1bSJed Brown { 177*071b71afSMatthew G. Knepley PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = {linear}; 178*071b71afSMatthew G. Knepley DM dmProj, dmaProj; 179c4762a1bSJed Brown PetscFE fe; 180*071b71afSMatthew G. Knepley KSP ksp; 181*071b71afSMatthew G. Knepley Mat Interp, mass, mass2; 182*071b71afSMatthew G. Knepley Vec u, ua, scaling, rhs, uproj; 183c4762a1bSJed Brown PetscReal error; 184*071b71afSMatthew G. Knepley PetscBool simplex; 185c4762a1bSJed Brown PetscInt dim; 186c4762a1bSJed Brown PetscErrorCode ierr; 187c4762a1bSJed Brown 188c4762a1bSJed Brown PetscFunctionBeginUser; 189c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 190*071b71afSMatthew G. Knepley ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 191c4762a1bSJed Brown 192*071b71afSMatthew G. Knepley ierr = DMClone(dm, &dmProj);CHKERRQ(ierr); 193*071b71afSMatthew G. Knepley ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe);CHKERRQ(ierr); 194*071b71afSMatthew G. Knepley ierr = DMSetField(dmProj, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 195*071b71afSMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 196*071b71afSMatthew G. Knepley ierr = DMCreateDS(dmProj);CHKERRQ(ierr); 197*071b71afSMatthew G. Knepley 198*071b71afSMatthew G. Knepley ierr = DMClone(dma, &dmaProj);CHKERRQ(ierr); 199*071b71afSMatthew G. Knepley ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe);CHKERRQ(ierr); 200*071b71afSMatthew G. Knepley ierr = DMSetField(dmaProj, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 201*071b71afSMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 202*071b71afSMatthew G. Knepley ierr = DMCreateDS(dmaProj);CHKERRQ(ierr); 203*071b71afSMatthew G. Knepley 204*071b71afSMatthew G. Knepley ierr = DMGetGlobalVector(dmProj, &u);CHKERRQ(ierr); 205*071b71afSMatthew G. Knepley ierr = DMGetGlobalVector(dmaProj, &ua);CHKERRQ(ierr); 206*071b71afSMatthew G. Knepley ierr = DMGetGlobalVector(dmaProj, &rhs);CHKERRQ(ierr); 207*071b71afSMatthew G. Knepley ierr = DMGetGlobalVector(dmaProj, &uproj);CHKERRQ(ierr); 208*071b71afSMatthew G. Knepley 209*071b71afSMatthew G. Knepley // Interpolate onto original mesh using dual basis 210*071b71afSMatthew G. Knepley ierr = DMProjectFunction(dmProj, 0.0, funcs, NULL, INSERT_VALUES, u);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "Original");CHKERRQ(ierr); 212c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-orig_vec_view");CHKERRQ(ierr); 213*071b71afSMatthew G. Knepley ierr = DMComputeL2Diff(dmProj, 0.0, funcs, NULL, u, &error);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Original L2 Error: %g\n", (double) error);CHKERRQ(ierr); 215*071b71afSMatthew G. Knepley // Interpolate onto NEW mesh using dual basis 216*071b71afSMatthew G. Knepley ierr = DMProjectFunction(dmaProj, 0.0, funcs, NULL, INSERT_VALUES, ua);CHKERRQ(ierr); 217*071b71afSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) ua, "Adapted");CHKERRQ(ierr); 218*071b71afSMatthew G. Knepley ierr = VecViewFromOptions(ua, NULL, "-adapt_vec_view");CHKERRQ(ierr); 219*071b71afSMatthew G. Knepley ierr = DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, ua, &error);CHKERRQ(ierr); 220*071b71afSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "Adapted L2 Error: %g\n", (double) error);CHKERRQ(ierr); 221*071b71afSMatthew G. Knepley // Interpolate between meshes using interpolation matrix 222*071b71afSMatthew G. Knepley ierr = DMCreateInterpolation(dmProj, dmaProj, &Interp, &scaling);CHKERRQ(ierr); 223c4762a1bSJed Brown ierr = MatInterpolate(Interp, u, ua);CHKERRQ(ierr); 224c4762a1bSJed Brown ierr = MatDestroy(&Interp);CHKERRQ(ierr); 225c4762a1bSJed Brown ierr = VecDestroy(&scaling);CHKERRQ(ierr); 226c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) ua, "Interpolation");CHKERRQ(ierr); 227c4762a1bSJed Brown ierr = VecViewFromOptions(ua, NULL, "-interp_vec_view");CHKERRQ(ierr); 228*071b71afSMatthew G. Knepley ierr = DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, ua, &error);CHKERRQ(ierr); 229c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Interpolated L2 Error: %g\n", (double) error);CHKERRQ(ierr); 230*071b71afSMatthew G. Knepley // L2 projection 231*071b71afSMatthew G. Knepley ierr = DMCreateMassMatrix(dmaProj, dmaProj, &mass);CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = MatViewFromOptions(mass, NULL, "-mass_mat_view");CHKERRQ(ierr); 233c4762a1bSJed Brown ierr = KSPCreate(PETSC_COMM_WORLD, &ksp);CHKERRQ(ierr); 234c4762a1bSJed Brown ierr = KSPSetOperators(ksp, mass, mass);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 236*071b71afSMatthew G. Knepley // Compute rhs as M f, could also direclty project the analytic function but we might not have it 237*071b71afSMatthew G. Knepley ierr = DMCreateMassMatrix(dmProj, dmaProj, &mass2);CHKERRQ(ierr); 238*071b71afSMatthew G. Knepley ierr = MatMult(mass2, u, rhs);CHKERRQ(ierr); 239*071b71afSMatthew G. Knepley ierr = MatDestroy(&mass2);CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = KSPSolve(ksp, rhs, uproj);CHKERRQ(ierr); 241*071b71afSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) uproj, "L_2 Projection");CHKERRQ(ierr); 242c4762a1bSJed Brown ierr = VecViewFromOptions(uproj, NULL, "-proj_vec_view");CHKERRQ(ierr); 243*071b71afSMatthew G. Knepley ierr = DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, uproj, &error);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double) error);CHKERRQ(ierr); 245*071b71afSMatthew G. Knepley ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 246c4762a1bSJed Brown ierr = MatDestroy(&mass);CHKERRQ(ierr); 247*071b71afSMatthew G. Knepley ierr = DMRestoreGlobalVector(dmProj, &u);CHKERRQ(ierr); 248*071b71afSMatthew G. Knepley ierr = DMRestoreGlobalVector(dmaProj, &ua);CHKERRQ(ierr); 249*071b71afSMatthew G. Knepley ierr = DMRestoreGlobalVector(dmaProj, &rhs);CHKERRQ(ierr); 250*071b71afSMatthew G. Knepley ierr = DMRestoreGlobalVector(dmaProj, &uproj);CHKERRQ(ierr); 251*071b71afSMatthew G. Knepley ierr = DMDestroy(&dmProj);CHKERRQ(ierr); 252*071b71afSMatthew G. Knepley ierr = DMDestroy(&dmaProj);CHKERRQ(ierr); 253c4762a1bSJed Brown PetscFunctionReturn(0); 254c4762a1bSJed Brown } 255c4762a1bSJed Brown 256c4762a1bSJed Brown int main (int argc, char * argv[]) { 257*071b71afSMatthew 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; 262*071b71afSMatthew G. Knepley PetscInt r; 263c4762a1bSJed Brown PetscErrorCode ierr; 264c4762a1bSJed Brown 265c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 266c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 267c4762a1bSJed Brown ierr = ProcessOptions(comm, &user);CHKERRQ(ierr); 268*071b71afSMatthew G. Knepley ierr = CreateMesh(comm, &dm);CHKERRQ(ierr); 269c4762a1bSJed Brown 270*071b71afSMatthew G. Knepley odm = dm; 271*071b71afSMatthew G. Knepley ierr = DMPlexDistributeOverlap(odm, 1, NULL, &dm);CHKERRQ(ierr); 272*071b71afSMatthew G. Knepley if (!dm) {dm = odm;} 273c4762a1bSJed Brown else {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 274*071b71afSMatthew G. Knepley 275*071b71afSMatthew G. Knepley for (r = 0; r < user.Nr; ++r) { 276*071b71afSMatthew G. Knepley DMLabel label; 277*071b71afSMatthew G. Knepley 278*071b71afSMatthew G. Knepley ierr = ComputeMetric(dm, &user, &metric);CHKERRQ(ierr); 279*071b71afSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 280*071b71afSMatthew G. Knepley ierr = DMAdaptMetric(dm, metric, label, &dma);CHKERRQ(ierr); 281*071b71afSMatthew G. Knepley ierr = VecDestroy(&metric);CHKERRQ(ierr); 282c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) dma, "DMadapt");CHKERRQ(ierr); 283c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject) dma, "adapt_");CHKERRQ(ierr); 284c4762a1bSJed Brown ierr = DMViewFromOptions(dma, NULL, "-dm_view");CHKERRQ(ierr); 285*071b71afSMatthew G. Knepley if (user.doL2) {ierr = TestL2Projection(dm, dma, &user);CHKERRQ(ierr);} 286*071b71afSMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 287*071b71afSMatthew G. Knepley dm = dma; 288*071b71afSMatthew G. Knepley } 289*071b71afSMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) dm, "final_");CHKERRQ(ierr); 290*071b71afSMatthew G. Knepley ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr); 291*071b71afSMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 292*071b71afSMatthew G. Knepley ierr = PetscFinalize(); 293*071b71afSMatthew G. Knepley return ierr; 294c4762a1bSJed Brown } 295c4762a1bSJed Brown 296c4762a1bSJed Brown /*TEST 297c4762a1bSJed Brown 2988e3a2eefSMatthew G. Knepley build: 2998e3a2eefSMatthew G. Knepley requires: pragmatic 3008e3a2eefSMatthew G. Knepley 301*071b71afSMatthew G. Knepley testset: 302*071b71afSMatthew G. Knepley args: -dm_plex_box_faces 4,4,4 -met 2 -init_dm_view -adapt_dm_view 303*071b71afSMatthew G. Knepley 304c4762a1bSJed Brown test: 305c4762a1bSJed Brown suffix: 0 306*071b71afSMatthew G. Knepley args: -dm_plex_separate_marker 0 307c4762a1bSJed Brown test: 308c4762a1bSJed Brown suffix: 1 309*071b71afSMatthew G. Knepley args: -dm_plex_separate_marker 1 310c4762a1bSJed Brown test: 311c4762a1bSJed Brown suffix: 2 312*071b71afSMatthew G. Knepley args: -dm_plex_dim 3 313c4762a1bSJed Brown test: 314c4762a1bSJed Brown suffix: 3 315*071b71afSMatthew G. Knepley args: -dm_plex_dim 3 316*071b71afSMatthew G. Knepley 317*071b71afSMatthew G. Knepley # Pragmatic hangs for simple partitioner 318*071b71afSMatthew G. Knepley testset: 319*071b71afSMatthew G. Knepley requires: parmetis 320*071b71afSMatthew G. Knepley args: -dm_plex_box_faces 2,2 -dm_distribute -petscpartitioner_type parmetis -met 2 -init_dm_view -adapt_dm_view 321*071b71afSMatthew G. Knepley 322c4762a1bSJed Brown test: 323c4762a1bSJed Brown suffix: 4 324c4762a1bSJed Brown nsize: 2 325c4762a1bSJed Brown test: 326c4762a1bSJed Brown suffix: 5 327c4762a1bSJed Brown nsize: 4 328*071b71afSMatthew G. Knepley 329c4762a1bSJed Brown test: 330c4762a1bSJed Brown suffix: 6 331c4762a1bSJed Brown nsize: 2 332*071b71afSMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_box_faces 9,9,9 -dm_distribute -petscpartitioner_type parmetis \ 333*071b71afSMatthew G. Knepley -met 0 -hmin 0.01 -hmax 0.03 -init_dm_view -adapt_dm_view 334c4762a1bSJed Brown test: 335c4762a1bSJed Brown suffix: 7 3363e32e87aSJoe Wallwork nsize: 2 337*071b71afSMatthew G. Knepley args: -dm_plex_box_faces 19,19 -dm_distribute -petscpartitioner_type parmetis \ 338*071b71afSMatthew G. Knepley -met 2 -hmax 0.5 -hmin 0.001 -init_dm_view -adapt_dm_view 339c4762a1bSJed Brown test: 340c4762a1bSJed Brown suffix: proj_0 341*071b71afSMatthew G. Knepley args: -dm_plex_box_faces 2,2 -dm_plex_hash_location -init_dm_view -adapt_dm_view -do_L2 \ 342*071b71afSMatthew G. Knepley -petscspace_degree 1 -petscfe_default_quadrature_order 1 -pc_type lu 343c4762a1bSJed Brown test: 344c4762a1bSJed Brown suffix: proj_1 345*071b71afSMatthew G. Knepley args: -dm_plex_box_faces 4,4 -dm_plex_hash_location -init_dm_view -adapt_dm_view -do_L2 \ 346*071b71afSMatthew G. Knepley -petscspace_degree 1 -petscfe_default_quadrature_order 1 -pc_type lu 347*071b71afSMatthew G. Knepley 348*071b71afSMatthew G. Knepley test: 349*071b71afSMatthew G. Knepley suffix: sensor 350*071b71afSMatthew G. Knepley args: -dm_plex_box_faces 9,9 -met 3 -init_dm_view -adapt_dm_view \ 351*071b71afSMatthew 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 \ 352*071b71afSMatthew G. Knepley -dm_plex_metric_target_complexity 10000.0 353c4762a1bSJed Brown 354c4762a1bSJed Brown TEST*/ 355