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 45c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Meshing Adaptation Options", "DMPLEX");CHKERRQ(ierr); 46071b71afSMatthew 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 56071b71afSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 57c4762a1bSJed Brown { 58c4762a1bSJed Brown PetscErrorCode ierr; 59c4762a1bSJed Brown 60c4762a1bSJed Brown PetscFunctionBegin; 61071b71afSMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 62071b71afSMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 63071b71afSMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 64071b71afSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *dm, "DMinit");CHKERRQ(ierr); 65071b71afSMatthew G. Knepley ierr = DMViewFromOptions(*dm, NULL, "-init_dm_view");CHKERRQ(ierr); 66071b71afSMatthew G. Knepley PetscFunctionReturn(0); 67c4762a1bSJed Brown } 68071b71afSMatthew G. Knepley 69071b71afSMatthew G. Knepley static PetscErrorCode ComputeMetricSensor(DM dm, AppCtx *user, Vec *metric) 70c4762a1bSJed Brown { 71071b71afSMatthew G. Knepley PetscSimplePointFunc funcs[1] = {sensor}; 72071b71afSMatthew G. Knepley DM dmSensor, dmGrad, dmHess; 73071b71afSMatthew G. Knepley PetscFE fe; 74071b71afSMatthew G. Knepley Vec f, g, H; 75071b71afSMatthew G. Knepley PetscBool simplex; 76071b71afSMatthew G. Knepley PetscInt dim; 77071b71afSMatthew G. Knepley PetscErrorCode ierr; 78c4762a1bSJed Brown 79071b71afSMatthew G. Knepley PetscFunctionBegin; 80071b71afSMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 81071b71afSMatthew G. Knepley ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 82071b71afSMatthew G. Knepley 83071b71afSMatthew G. Knepley ierr = DMClone(dm, &dmSensor);CHKERRQ(ierr); 84071b71afSMatthew G. Knepley ierr = PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, simplex, 1, -1, &fe);CHKERRQ(ierr); 85071b71afSMatthew G. Knepley ierr = DMSetField(dmSensor, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 86071b71afSMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 87071b71afSMatthew G. Knepley ierr = DMCreateDS(dmSensor);CHKERRQ(ierr); 88071b71afSMatthew G. Knepley ierr = DMCreateLocalVector(dmSensor, &f);CHKERRQ(ierr); 89071b71afSMatthew G. Knepley ierr = DMProjectFunctionLocal(dmSensor, 0., funcs, NULL, INSERT_VALUES, f);CHKERRQ(ierr); 90071b71afSMatthew G. Knepley ierr = VecViewFromOptions(f, NULL, "-sensor_view");CHKERRQ(ierr); 91071b71afSMatthew G. Knepley 92071b71afSMatthew G. Knepley // Recover the gradient of the sensor function 93071b71afSMatthew G. Knepley ierr = DMClone(dm, &dmGrad);CHKERRQ(ierr); 94071b71afSMatthew G. Knepley ierr = PetscFECreateLagrange(PETSC_COMM_SELF, dim, dim, simplex, 1, -1, &fe);CHKERRQ(ierr); 95071b71afSMatthew G. Knepley ierr = DMSetField(dmGrad, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 96071b71afSMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 97071b71afSMatthew G. Knepley ierr = DMCreateDS(dmGrad);CHKERRQ(ierr); 98071b71afSMatthew G. Knepley ierr = DMCreateLocalVector(dmGrad, &g);CHKERRQ(ierr); 99071b71afSMatthew G. Knepley ierr = DMPlexComputeGradientClementInterpolant(dmSensor, f, g);CHKERRQ(ierr); 100071b71afSMatthew G. Knepley ierr = VecDestroy(&f);CHKERRQ(ierr); 101071b71afSMatthew G. Knepley ierr = VecViewFromOptions(g, NULL, "-gradient_view");CHKERRQ(ierr); 102071b71afSMatthew G. Knepley 103071b71afSMatthew G. Knepley // Recover the Hessian of the sensor function 104071b71afSMatthew G. Knepley ierr = DMClone(dm, &dmHess);CHKERRQ(ierr); 105071b71afSMatthew G. Knepley ierr = DMPlexMetricCreate(dmHess, 0, &H);CHKERRQ(ierr); 106071b71afSMatthew G. Knepley ierr = DMPlexComputeGradientClementInterpolant(dmGrad, g, H);CHKERRQ(ierr); 107071b71afSMatthew G. Knepley ierr = VecDestroy(&g);CHKERRQ(ierr); 108071b71afSMatthew G. Knepley ierr = VecViewFromOptions(H, NULL, "-hessian_view");CHKERRQ(ierr); 109071b71afSMatthew G. Knepley 110071b71afSMatthew G. Knepley // Obtain a metric by Lp normalization 111071b71afSMatthew G. Knepley ierr = DMPlexMetricNormalize(dmHess, H, PETSC_TRUE, PETSC_TRUE, metric);CHKERRQ(ierr); 112071b71afSMatthew G. Knepley ierr = VecDestroy(&H);CHKERRQ(ierr); 113071b71afSMatthew G. Knepley ierr = DMDestroy(&dmHess);CHKERRQ(ierr); 114071b71afSMatthew G. Knepley ierr = DMDestroy(&dmGrad);CHKERRQ(ierr); 115071b71afSMatthew 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 { 121071b71afSMatthew G. Knepley PetscReal lambda = 1/(user->hmax*user->hmax); 122c4762a1bSJed Brown PetscErrorCode ierr; 123c4762a1bSJed Brown 124c4762a1bSJed Brown PetscFunctionBeginUser; 125071b71afSMatthew G. Knepley if (user->metOpt == 0) { 1263e32e87aSJoe Wallwork /* Specify a uniform, isotropic metric */ 1273e32e87aSJoe Wallwork ierr = DMPlexMetricCreateUniform(dm, 0, lambda, metric);CHKERRQ(ierr); 128071b71afSMatthew G. Knepley } else if (user->metOpt == 3) { 129071b71afSMatthew G. Knepley ierr = ComputeMetricSensor(dm, user, metric);CHKERRQ(ierr); 130071b71afSMatthew 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*fdc00aefSJoe Wallwork PetscInt dim, i, j, vStart, vEnd, v; 1373e32e87aSJoe Wallwork 138*fdc00aefSJoe Wallwork ierr = DMPlexMetricCreate(dm, 0, metric);CHKERRQ(ierr); 139*fdc00aefSJoe Wallwork ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = VecGetArray(*metric, &met);CHKERRQ(ierr); 144071b71afSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 145071b71afSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 146071b71afSMatthew G. Knepley PetscScalar *vcoords; 147c4762a1bSJed Brown PetscScalar *pmet; 148c4762a1bSJed Brown 149071b71afSMatthew G. Knepley ierr = DMPlexPointLocalRead(cdm, v, coords, &vcoords);CHKERRQ(ierr); 150c4762a1bSJed Brown switch (user->metOpt) { 151c4762a1bSJed Brown case 1: 152071b71afSMatthew G. Knepley h = user->hmax - (user->hmax-user->hmin)*PetscRealPart(vcoords[0]); 153c4762a1bSJed Brown break; 154c4762a1bSJed Brown case 2: 155071b71afSMatthew G. Knepley h = user->hmax*PetscAbsReal(((PetscReal) 1.0)-PetscExpReal(-PetscAbsScalar(vcoords[0]-(PetscReal)0.5))) + user->hmin; 156c4762a1bSJed Brown break; 157c4762a1bSJed Brown default: 158*fdc00aefSJoe Wallwork SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "metOpt = 0, 1, 2 or 3, cannot be %d", user->metOpt); 159c4762a1bSJed Brown } 160071b71afSMatthew G. Knepley ierr = DMPlexPointLocalRef(dm, v, met, &pmet);CHKERRQ(ierr); 161*fdc00aefSJoe Wallwork for (i = 0; i < dim; ++i) { 162*fdc00aefSJoe Wallwork for (j = 0; j < dim; ++j) { 163*fdc00aefSJoe Wallwork if (i == j) { 164*fdc00aefSJoe Wallwork if (i == 0) pmet[i*dim+j] = 1/(h*h); 165*fdc00aefSJoe Wallwork else pmet[i*dim+j] = lambda; 166*fdc00aefSJoe Wallwork } else pmet[i*dim+j] = 0.0; 167*fdc00aefSJoe Wallwork } 168*fdc00aefSJoe Wallwork } 169c4762a1bSJed Brown } 170c4762a1bSJed Brown ierr = VecRestoreArray(*metric, &met);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 1723e32e87aSJoe Wallwork } 173c4762a1bSJed Brown PetscFunctionReturn(0); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown 176c4762a1bSJed Brown static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 177c4762a1bSJed Brown { 178c4762a1bSJed Brown u[0] = x[0] + x[1]; 179c4762a1bSJed Brown return 0; 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown static PetscErrorCode TestL2Projection(DM dm, DM dma, AppCtx *user) 183c4762a1bSJed Brown { 184071b71afSMatthew G. Knepley PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *) = {linear}; 185071b71afSMatthew G. Knepley DM dmProj, dmaProj; 186c4762a1bSJed Brown PetscFE fe; 187071b71afSMatthew G. Knepley KSP ksp; 188071b71afSMatthew G. Knepley Mat Interp, mass, mass2; 189071b71afSMatthew G. Knepley Vec u, ua, scaling, rhs, uproj; 190c4762a1bSJed Brown PetscReal error; 191071b71afSMatthew G. Knepley PetscBool simplex; 192c4762a1bSJed Brown PetscInt dim; 193c4762a1bSJed Brown PetscErrorCode ierr; 194c4762a1bSJed Brown 195c4762a1bSJed Brown PetscFunctionBeginUser; 196c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 197071b71afSMatthew G. Knepley ierr = DMPlexIsSimplex(dm, &simplex);CHKERRQ(ierr); 198c4762a1bSJed Brown 199071b71afSMatthew G. Knepley ierr = DMClone(dm, &dmProj);CHKERRQ(ierr); 200071b71afSMatthew G. Knepley ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe);CHKERRQ(ierr); 201071b71afSMatthew G. Knepley ierr = DMSetField(dmProj, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 202071b71afSMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 203071b71afSMatthew G. Knepley ierr = DMCreateDS(dmProj);CHKERRQ(ierr); 204071b71afSMatthew G. Knepley 205071b71afSMatthew G. Knepley ierr = DMClone(dma, &dmaProj);CHKERRQ(ierr); 206071b71afSMatthew G. Knepley ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, 1, simplex, NULL, -1, &fe);CHKERRQ(ierr); 207071b71afSMatthew G. Knepley ierr = DMSetField(dmaProj, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 208071b71afSMatthew G. Knepley ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 209071b71afSMatthew G. Knepley ierr = DMCreateDS(dmaProj);CHKERRQ(ierr); 210071b71afSMatthew G. Knepley 211071b71afSMatthew G. Knepley ierr = DMGetGlobalVector(dmProj, &u);CHKERRQ(ierr); 212071b71afSMatthew G. Knepley ierr = DMGetGlobalVector(dmaProj, &ua);CHKERRQ(ierr); 213071b71afSMatthew G. Knepley ierr = DMGetGlobalVector(dmaProj, &rhs);CHKERRQ(ierr); 214071b71afSMatthew G. Knepley ierr = DMGetGlobalVector(dmaProj, &uproj);CHKERRQ(ierr); 215071b71afSMatthew G. Knepley 216071b71afSMatthew G. Knepley // Interpolate onto original mesh using dual basis 217071b71afSMatthew G. Knepley ierr = DMProjectFunction(dmProj, 0.0, funcs, NULL, INSERT_VALUES, u);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "Original");CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-orig_vec_view");CHKERRQ(ierr); 220071b71afSMatthew G. Knepley ierr = DMComputeL2Diff(dmProj, 0.0, funcs, NULL, u, &error);CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Original L2 Error: %g\n", (double) error);CHKERRQ(ierr); 222071b71afSMatthew G. Knepley // Interpolate onto NEW mesh using dual basis 223071b71afSMatthew G. Knepley ierr = DMProjectFunction(dmaProj, 0.0, funcs, NULL, INSERT_VALUES, ua);CHKERRQ(ierr); 224071b71afSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) ua, "Adapted");CHKERRQ(ierr); 225071b71afSMatthew G. Knepley ierr = VecViewFromOptions(ua, NULL, "-adapt_vec_view");CHKERRQ(ierr); 226071b71afSMatthew G. Knepley ierr = DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, ua, &error);CHKERRQ(ierr); 227071b71afSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "Adapted L2 Error: %g\n", (double) error);CHKERRQ(ierr); 228071b71afSMatthew G. Knepley // Interpolate between meshes using interpolation matrix 229071b71afSMatthew G. Knepley ierr = DMCreateInterpolation(dmProj, dmaProj, &Interp, &scaling);CHKERRQ(ierr); 230c4762a1bSJed Brown ierr = MatInterpolate(Interp, u, ua);CHKERRQ(ierr); 231c4762a1bSJed Brown ierr = MatDestroy(&Interp);CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = VecDestroy(&scaling);CHKERRQ(ierr); 233c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) ua, "Interpolation");CHKERRQ(ierr); 234c4762a1bSJed Brown ierr = VecViewFromOptions(ua, NULL, "-interp_vec_view");CHKERRQ(ierr); 235071b71afSMatthew G. Knepley ierr = DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, ua, &error);CHKERRQ(ierr); 236c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Interpolated L2 Error: %g\n", (double) error);CHKERRQ(ierr); 237071b71afSMatthew G. Knepley // L2 projection 238071b71afSMatthew G. Knepley ierr = DMCreateMassMatrix(dmaProj, dmaProj, &mass);CHKERRQ(ierr); 239c4762a1bSJed Brown ierr = MatViewFromOptions(mass, NULL, "-mass_mat_view");CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = KSPCreate(PETSC_COMM_WORLD, &ksp);CHKERRQ(ierr); 241c4762a1bSJed Brown ierr = KSPSetOperators(ksp, mass, mass);CHKERRQ(ierr); 242c4762a1bSJed Brown ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 243071b71afSMatthew G. Knepley // Compute rhs as M f, could also direclty project the analytic function but we might not have it 244071b71afSMatthew G. Knepley ierr = DMCreateMassMatrix(dmProj, dmaProj, &mass2);CHKERRQ(ierr); 245071b71afSMatthew G. Knepley ierr = MatMult(mass2, u, rhs);CHKERRQ(ierr); 246071b71afSMatthew G. Knepley ierr = MatDestroy(&mass2);CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = KSPSolve(ksp, rhs, uproj);CHKERRQ(ierr); 248071b71afSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) uproj, "L_2 Projection");CHKERRQ(ierr); 249c4762a1bSJed Brown ierr = VecViewFromOptions(uproj, NULL, "-proj_vec_view");CHKERRQ(ierr); 250071b71afSMatthew G. Knepley ierr = DMComputeL2Diff(dmaProj, 0.0, funcs, NULL, uproj, &error);CHKERRQ(ierr); 251c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double) error);CHKERRQ(ierr); 252071b71afSMatthew G. Knepley ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 253c4762a1bSJed Brown ierr = MatDestroy(&mass);CHKERRQ(ierr); 254071b71afSMatthew G. Knepley ierr = DMRestoreGlobalVector(dmProj, &u);CHKERRQ(ierr); 255071b71afSMatthew G. Knepley ierr = DMRestoreGlobalVector(dmaProj, &ua);CHKERRQ(ierr); 256071b71afSMatthew G. Knepley ierr = DMRestoreGlobalVector(dmaProj, &rhs);CHKERRQ(ierr); 257071b71afSMatthew G. Knepley ierr = DMRestoreGlobalVector(dmaProj, &uproj);CHKERRQ(ierr); 258071b71afSMatthew G. Knepley ierr = DMDestroy(&dmProj);CHKERRQ(ierr); 259071b71afSMatthew G. Knepley ierr = DMDestroy(&dmaProj);CHKERRQ(ierr); 260c4762a1bSJed Brown PetscFunctionReturn(0); 261c4762a1bSJed Brown } 262c4762a1bSJed Brown 263c4762a1bSJed Brown int main (int argc, char * argv[]) { 264071b71afSMatthew G. Knepley DM dm; 265c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 266c4762a1bSJed Brown MPI_Comm comm; 267c4762a1bSJed Brown DM dma, odm; 268c4762a1bSJed Brown Vec metric; 269071b71afSMatthew G. Knepley PetscInt r; 270c4762a1bSJed Brown PetscErrorCode ierr; 271c4762a1bSJed Brown 272c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 273c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 274c4762a1bSJed Brown ierr = ProcessOptions(comm, &user);CHKERRQ(ierr); 275071b71afSMatthew G. Knepley ierr = CreateMesh(comm, &dm);CHKERRQ(ierr); 276c4762a1bSJed Brown 277071b71afSMatthew G. Knepley odm = dm; 278071b71afSMatthew G. Knepley ierr = DMPlexDistributeOverlap(odm, 1, NULL, &dm);CHKERRQ(ierr); 279071b71afSMatthew G. Knepley if (!dm) {dm = odm;} 280c4762a1bSJed Brown else {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 281071b71afSMatthew G. Knepley 282071b71afSMatthew G. Knepley for (r = 0; r < user.Nr; ++r) { 283071b71afSMatthew G. Knepley DMLabel label; 284071b71afSMatthew G. Knepley 285071b71afSMatthew G. Knepley ierr = ComputeMetric(dm, &user, &metric);CHKERRQ(ierr); 286071b71afSMatthew G. Knepley ierr = DMGetLabel(dm, "marker", &label);CHKERRQ(ierr); 287e9bc5bffSJoe Wallwork ierr = DMAdaptMetric(dm, metric, label, NULL, &dma);CHKERRQ(ierr); 288071b71afSMatthew G. Knepley ierr = VecDestroy(&metric);CHKERRQ(ierr); 289c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) dma, "DMadapt");CHKERRQ(ierr); 290c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject) dma, "adapt_");CHKERRQ(ierr); 291c4762a1bSJed Brown ierr = DMViewFromOptions(dma, NULL, "-dm_view");CHKERRQ(ierr); 292071b71afSMatthew G. Knepley if (user.doL2) {ierr = TestL2Projection(dm, dma, &user);CHKERRQ(ierr);} 293071b71afSMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 294071b71afSMatthew G. Knepley dm = dma; 295071b71afSMatthew G. Knepley } 296071b71afSMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) dm, "final_");CHKERRQ(ierr); 297071b71afSMatthew G. Knepley ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr); 298071b71afSMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 299071b71afSMatthew G. Knepley ierr = PetscFinalize(); 300071b71afSMatthew G. Knepley return ierr; 301c4762a1bSJed Brown } 302c4762a1bSJed Brown 303c4762a1bSJed Brown /*TEST 304c4762a1bSJed Brown 3058e3a2eefSMatthew G. Knepley build: 3068e3a2eefSMatthew G. Knepley requires: pragmatic 3078e3a2eefSMatthew G. Knepley 308071b71afSMatthew G. Knepley testset: 3099e3182e6SJoe Wallwork args: -dm_plex_box_faces 4,4,4 -dm_adaptor pragmatic -met 2 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic 310071b71afSMatthew G. Knepley 311c4762a1bSJed Brown test: 312c4762a1bSJed Brown suffix: 0 313071b71afSMatthew G. Knepley args: -dm_plex_separate_marker 0 314c4762a1bSJed Brown test: 315c4762a1bSJed Brown suffix: 1 316071b71afSMatthew G. Knepley args: -dm_plex_separate_marker 1 317c4762a1bSJed Brown test: 318c4762a1bSJed Brown suffix: 2 319071b71afSMatthew G. Knepley args: -dm_plex_dim 3 320c4762a1bSJed Brown test: 321c4762a1bSJed Brown suffix: 3 322071b71afSMatthew G. Knepley args: -dm_plex_dim 3 323071b71afSMatthew G. Knepley 324071b71afSMatthew G. Knepley # Pragmatic hangs for simple partitioner 325071b71afSMatthew G. Knepley testset: 326071b71afSMatthew G. Knepley requires: parmetis 3279e3182e6SJoe Wallwork args: -dm_plex_box_faces 2,2 -dm_adaptor pragmatic -dm_distribute -petscpartitioner_type parmetis -met 2 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic 328071b71afSMatthew G. Knepley 329c4762a1bSJed Brown test: 330c4762a1bSJed Brown suffix: 4 331c4762a1bSJed Brown nsize: 2 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown suffix: 5 334c4762a1bSJed Brown nsize: 4 335071b71afSMatthew G. Knepley 336c4762a1bSJed Brown test: 3379880b877SBarry Smith requires: parmetis 338c4762a1bSJed Brown suffix: 6 339c4762a1bSJed Brown nsize: 2 340cbddfcd8SJoe Wallwork args: -dm_plex_dim 3 -dm_plex_box_faces 9,9,9 -dm_adaptor pragmatic -dm_distribute -petscpartitioner_type parmetis \ 3419e3182e6SJoe Wallwork -met 0 -hmin 0.01 -hmax 0.03 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic 342c4762a1bSJed Brown test: 3439880b877SBarry Smith requires: parmetis 344c4762a1bSJed Brown suffix: 7 3453e32e87aSJoe Wallwork nsize: 2 346cbddfcd8SJoe Wallwork args: -dm_plex_box_faces 19,19 -dm_adaptor pragmatic -dm_distribute -petscpartitioner_type parmetis \ 3479e3182e6SJoe Wallwork -met 2 -hmax 0.5 -hmin 0.001 -init_dm_view -adapt_dm_view -dm_adaptor pragmatic 348c4762a1bSJed Brown test: 349c4762a1bSJed Brown suffix: proj_0 350cbddfcd8SJoe Wallwork args: -dm_plex_box_faces 2,2 -dm_plex_hash_location -dm_adaptor pragmatic -init_dm_view -adapt_dm_view -do_L2 \ 3519e3182e6SJoe Wallwork -petscspace_degree 1 -petscfe_default_quadrature_order 1 -pc_type lu -dm_adaptor pragmatic 352c4762a1bSJed Brown test: 353c4762a1bSJed Brown suffix: proj_1 354cbddfcd8SJoe Wallwork args: -dm_plex_box_faces 4,4 -dm_plex_hash_location -dm_adaptor pragmatic -init_dm_view -adapt_dm_view -do_L2 \ 3559e3182e6SJoe Wallwork -petscspace_degree 1 -petscfe_default_quadrature_order 1 -pc_type lu -dm_adaptor pragmatic 356071b71afSMatthew G. Knepley 357071b71afSMatthew G. Knepley test: 358071b71afSMatthew G. Knepley suffix: sensor 359cbddfcd8SJoe Wallwork args: -dm_plex_box_faces 9,9 -met 3 -dm_adaptor pragmatic -init_dm_view -adapt_dm_view \ 360071b71afSMatthew 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 \ 3619e3182e6SJoe Wallwork -dm_plex_metric_target_complexity 10000.0 -dm_adaptor pragmatic 362c4762a1bSJed Brown 363c4762a1bSJed Brown TEST*/ 364