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 { 8c4762a1bSJed Brown DM dm; 9c4762a1bSJed Brown /* Definition of the test case (mesh and metric field) */ 10c4762a1bSJed Brown PetscInt dim; /* The topological mesh dimension */ 11c4762a1bSJed Brown char mshNam[PETSC_MAX_PATH_LEN]; /* Name of the mesh filename if any */ 12c4762a1bSJed Brown PetscInt nbrVerEdge; /* Number of vertices per edge if unit square/cube generated */ 13c4762a1bSJed Brown char bdLabel[PETSC_MAX_PATH_LEN]; /* Name of the label marking boundary facets */ 14c4762a1bSJed Brown PetscInt metOpt; /* Different choices of metric */ 15c4762a1bSJed Brown PetscReal hmax, hmin; /* Max and min sizes prescribed by the metric */ 16c4762a1bSJed Brown PetscBool doL2; /* Test L2 projection */ 17c4762a1bSJed Brown } AppCtx; 18c4762a1bSJed Brown 19c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 20c4762a1bSJed Brown { 21c4762a1bSJed Brown PetscErrorCode ierr; 22c4762a1bSJed Brown 23c4762a1bSJed Brown PetscFunctionBegin; 24c4762a1bSJed Brown options->dim = 2; 25c4762a1bSJed Brown ierr = PetscStrcpy(options->mshNam, "");CHKERRQ(ierr); 26c4762a1bSJed Brown options->nbrVerEdge = 5; 27c4762a1bSJed Brown ierr = PetscStrcpy(options->bdLabel, "");CHKERRQ(ierr); 28c4762a1bSJed Brown options->metOpt = 1; 29c4762a1bSJed Brown options->hmin = 0.05; 30c4762a1bSJed Brown options->hmax = 0.5; 31c4762a1bSJed Brown options->doL2 = PETSC_FALSE; 32c4762a1bSJed Brown 33c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "Meshing Adaptation Options", "DMPLEX");CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex19.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr); 358e3a2eefSMatthew G. Knepley ierr = PetscOptionsString("-msh", "Name of the mesh filename if any", "ex19.c", options->mshNam, options->mshNam, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr); 36c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-nbrVerEdge", "Number of vertices per edge if unit square/cube generated", "ex19.c", options->nbrVerEdge, &options->nbrVerEdge, NULL,0);CHKERRQ(ierr); 378e3a2eefSMatthew G. Knepley ierr = PetscOptionsString("-bdLabel", "Name of the label marking boundary facets", "ex19.c", options->bdLabel, options->bdLabel, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = PetscOptionsBoundedInt("-met", "Different choices of metric", "ex19.c", options->metOpt, &options->metOpt, NULL,0);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = PetscOptionsReal("-hmax", "Max size prescribed by the metric", "ex19.c", options->hmax, &options->hmax, NULL);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = PetscOptionsReal("-hmin", "Min size prescribed by the metric", "ex19.c", options->hmin, &options->hmin, NULL);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = PetscOptionsBool("-do_L2", "Test L2 projection", "ex19.c", options->doL2, &options->doL2, NULL);CHKERRQ(ierr); 421e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 43c4762a1bSJed Brown 44c4762a1bSJed Brown PetscFunctionReturn(0); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user) 48c4762a1bSJed Brown { 49c4762a1bSJed Brown PetscBool flag; 50c4762a1bSJed Brown PetscErrorCode ierr; 51c4762a1bSJed Brown 52c4762a1bSJed Brown PetscFunctionBegin; 53c4762a1bSJed Brown ierr = PetscStrcmp(user->mshNam, "", &flag);CHKERRQ(ierr); 54c4762a1bSJed Brown if (flag) { 55c4762a1bSJed Brown PetscInt faces[3]; 56c4762a1bSJed Brown faces[0] = user->nbrVerEdge-1; 57c4762a1bSJed Brown faces[1] = user->nbrVerEdge-1; 58c4762a1bSJed Brown faces[2] = user->nbrVerEdge-1; 59c4762a1bSJed Brown 60c4762a1bSJed Brown ierr = DMPlexCreateBoxMesh(comm, user->dim, PETSC_TRUE, faces, NULL, NULL, NULL, PETSC_TRUE, &user->dm);CHKERRQ(ierr); 61c4762a1bSJed Brown } else { 62c4762a1bSJed Brown ierr = DMPlexCreateFromFile(comm, user->mshNam, PETSC_TRUE, &user->dm);CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = DMGetDimension(user->dm, &user->dim);CHKERRQ(ierr); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown { 66c4762a1bSJed Brown DM distributedMesh = NULL; 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Distribute mesh over processes */ 69c4762a1bSJed Brown ierr = DMPlexDistribute(user->dm, 0, NULL, &distributedMesh);CHKERRQ(ierr); 70c4762a1bSJed Brown if (distributedMesh) { 71c4762a1bSJed Brown ierr = DMDestroy(&user->dm);CHKERRQ(ierr); 72c4762a1bSJed Brown user->dm = distributedMesh; 73c4762a1bSJed Brown } 74c4762a1bSJed Brown } 75c4762a1bSJed Brown ierr = DMSetFromOptions(user->dm);CHKERRQ(ierr); 76c4762a1bSJed Brown PetscFunctionReturn(0); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79c4762a1bSJed Brown static PetscErrorCode ComputeMetric(DM dm, AppCtx *user, Vec *metric) 80c4762a1bSJed Brown { 81*3e32e87aSJoe Wallwork PetscReal lambda; 82c4762a1bSJed Brown PetscErrorCode ierr; 83c4762a1bSJed Brown 84c4762a1bSJed Brown PetscFunctionBeginUser; 85*3e32e87aSJoe Wallwork 86*3e32e87aSJoe Wallwork /* Specify a uniform, isotropic metric */ 87*3e32e87aSJoe Wallwork lambda = 1/(user->hmax*user->hmax); 88*3e32e87aSJoe Wallwork ierr = DMPlexMetricCreateUniform(dm, 0, lambda, metric);CHKERRQ(ierr); 89*3e32e87aSJoe Wallwork 90*3e32e87aSJoe Wallwork /* Modify the first diagonal entry, if requested */ 91*3e32e87aSJoe Wallwork if (user->metOpt != 0) { 92*3e32e87aSJoe Wallwork DM cdm; 93*3e32e87aSJoe Wallwork PetscSection csec; 94*3e32e87aSJoe Wallwork Vec coordinates; 95*3e32e87aSJoe Wallwork const PetscScalar *coords; 96*3e32e87aSJoe Wallwork PetscScalar *met; 97*3e32e87aSJoe Wallwork PetscReal h; 98*3e32e87aSJoe Wallwork PetscInt pStart, pEnd, p; 99*3e32e87aSJoe Wallwork 100c4762a1bSJed Brown ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 101c4762a1bSJed Brown ierr = DMGetLocalSection(cdm, &csec);CHKERRQ(ierr); 102c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 103c4762a1bSJed Brown ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = VecGetArray(*metric, &met);CHKERRQ(ierr); 105*3e32e87aSJoe Wallwork ierr = PetscSectionGetChart(csec, &pStart, &pEnd);CHKERRQ(ierr); 106c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 107c4762a1bSJed Brown PetscScalar *pcoords; 108c4762a1bSJed Brown PetscScalar *pmet; 109c4762a1bSJed Brown 110c4762a1bSJed Brown ierr = DMPlexPointLocalRead(cdm, p, coords, &pcoords);CHKERRQ(ierr); 111c4762a1bSJed Brown switch (user->metOpt) { 112c4762a1bSJed Brown case 1: 113c4762a1bSJed Brown h = user->hmax - (user->hmax-user->hmin)*PetscRealPart(pcoords[0]); 114c4762a1bSJed Brown break; 115c4762a1bSJed Brown case 2: 116c4762a1bSJed Brown h = user->hmax*PetscAbsReal(((PetscReal) 1.0)-PetscExpReal(-PetscAbsScalar(pcoords[0]-(PetscReal)0.5))) + user->hmin; 117c4762a1bSJed Brown break; 118c4762a1bSJed Brown default: 119c4762a1bSJed Brown SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "metOpt = 0, 1 or 2, cannot be %d", user->metOpt); 120c4762a1bSJed Brown } 121*3e32e87aSJoe Wallwork lambda = 1/(h*h); 122*3e32e87aSJoe Wallwork ierr = DMPlexPointLocalRef(dm, p, met, &pmet);CHKERRQ(ierr); 123*3e32e87aSJoe Wallwork pmet[0] = lambda; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown ierr = VecRestoreArray(*metric, &met);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 127*3e32e87aSJoe Wallwork } 128c4762a1bSJed Brown PetscFunctionReturn(0); 129c4762a1bSJed Brown } 130c4762a1bSJed Brown 131c4762a1bSJed Brown static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 132c4762a1bSJed Brown { 133c4762a1bSJed Brown u[0] = x[0] + x[1]; 134c4762a1bSJed Brown return 0; 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137c4762a1bSJed Brown static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, 138c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 139c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 140c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 141c4762a1bSJed Brown { 142c4762a1bSJed Brown g0[0] = 1.0; 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145c4762a1bSJed Brown static PetscErrorCode TestL2Projection(DM dm, DM dma, AppCtx *user) 146c4762a1bSJed Brown { 147c4762a1bSJed Brown PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 148c4762a1bSJed Brown KSP ksp; 149c4762a1bSJed Brown PetscDS prob; 150c4762a1bSJed Brown PetscFE fe; 151c4762a1bSJed Brown Mat Interp, mass; 152c4762a1bSJed Brown Vec u, ua, scaling, ones, massLumped, rhs, uproj; 153c4762a1bSJed Brown PetscReal error; 154c4762a1bSJed Brown PetscInt dim; 155c4762a1bSJed Brown MPI_Comm comm; 156c4762a1bSJed Brown PetscErrorCode ierr; 157c4762a1bSJed Brown 158c4762a1bSJed Brown PetscFunctionBeginUser; 159c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, NULL, -1, &fe);CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = PetscDSSetDiscretization(prob, 0, (PetscObject) fe);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, NULL, -1, &fe);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = DMGetDS(dma, &prob);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = PetscDSSetDiscretization(prob, 0, (PetscObject) fe);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, identity, NULL, NULL, NULL);CHKERRQ(ierr); 170c4762a1bSJed Brown 171c4762a1bSJed Brown funcs[0] = linear; 172c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr); 173c4762a1bSJed Brown ierr = DMGetGlobalVector(dma, &ua);CHKERRQ(ierr); 174c4762a1bSJed Brown ierr = DMGetGlobalVector(dma, &ones);CHKERRQ(ierr); 175c4762a1bSJed Brown ierr = DMGetGlobalVector(dma, &massLumped);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = DMGetGlobalVector(dma, &rhs);CHKERRQ(ierr); 177c4762a1bSJed Brown ierr = DMGetGlobalVector(dma, &uproj);CHKERRQ(ierr); 178c4762a1bSJed Brown ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_VALUES, u);CHKERRQ(ierr); 179c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "Original");CHKERRQ(ierr); 180c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-orig_vec_view");CHKERRQ(ierr); 181c4762a1bSJed Brown ierr = DMComputeL2Diff(dm, 0.0, funcs, NULL, u, &error);CHKERRQ(ierr); 182c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Original L2 Error: %g\n", (double) error);CHKERRQ(ierr); 183c4762a1bSJed Brown ierr = DMCreateInterpolation(dm, dma, &Interp, &scaling);CHKERRQ(ierr); 184c4762a1bSJed Brown ierr = MatInterpolate(Interp, u, ua);CHKERRQ(ierr); 185c4762a1bSJed Brown ierr = MatDestroy(&Interp);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = VecDestroy(&scaling);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) ua, "Interpolation");CHKERRQ(ierr); 188c4762a1bSJed Brown ierr = VecViewFromOptions(ua, NULL, "-interp_vec_view");CHKERRQ(ierr); 189c4762a1bSJed Brown ierr = DMComputeL2Diff(dma, 0.0, funcs, NULL, ua, &error);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Interpolated L2 Error: %g\n", (double) error);CHKERRQ(ierr); 191c4762a1bSJed Brown 192c4762a1bSJed Brown ierr = VecSet(ones, 1.0);CHKERRQ(ierr); 1938e3a2eefSMatthew G. Knepley ierr = DMSNESComputeJacobianAction(dma, ua, ones, massLumped, user);CHKERRQ(ierr); 194c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) massLumped, "Lumped mass");CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = VecViewFromOptions(massLumped, NULL, "-mass_vec_view");CHKERRQ(ierr); 196c4762a1bSJed Brown ierr = DMCreateMassMatrix(dm, dma, &mass);CHKERRQ(ierr); 197c4762a1bSJed Brown ierr = MatMult(mass, u, rhs);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = MatDestroy(&mass);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = VecViewFromOptions(rhs, NULL, "-lumped_rhs_view");CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = VecPointwiseDivide(uproj, rhs, massLumped);CHKERRQ(ierr); 201c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) uproj, "Different Lumped Projection");CHKERRQ(ierr); 202c4762a1bSJed Brown ierr = VecViewFromOptions(uproj, NULL, "-lumped_rhs_vec_view");CHKERRQ(ierr); 203c4762a1bSJed Brown ierr = DMComputeL2Diff(dma, 0.0, funcs, NULL, uproj, &error);CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Lumped (rhs) L2 Error: %g\n", (double) error);CHKERRQ(ierr); 205c4762a1bSJed Brown 206c4762a1bSJed Brown ierr = DMCreateMatrix(dma, &mass);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = DMPlexSNESComputeJacobianFEM(dma, ua, mass, mass, user);CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = MatViewFromOptions(mass, NULL, "-mass_mat_view");CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = KSPCreate(PETSC_COMM_WORLD, &ksp);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = KSPSetOperators(ksp, mass, mass);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 212c4762a1bSJed Brown ierr = KSPSolve(ksp, rhs, uproj);CHKERRQ(ierr); 213c4762a1bSJed Brown ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) uproj, "Full Projection");CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = VecViewFromOptions(uproj, NULL, "-proj_vec_view");CHKERRQ(ierr); 216c4762a1bSJed Brown ierr = DMComputeL2Diff(dma, 0.0, funcs, NULL, uproj, &error);CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double) error);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = MatDestroy(&mass);CHKERRQ(ierr); 219c4762a1bSJed Brown 220c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dma, &ua);CHKERRQ(ierr); 222c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dma, &ones);CHKERRQ(ierr); 223c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dma, &massLumped);CHKERRQ(ierr); 224c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dma, &rhs);CHKERRQ(ierr); 225c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dma, &uproj);CHKERRQ(ierr); 226c4762a1bSJed Brown PetscFunctionReturn(0); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown 229c4762a1bSJed Brown int main (int argc, char * argv[]) { 230c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 231c4762a1bSJed Brown DMLabel bdLabel = NULL; 232c4762a1bSJed Brown MPI_Comm comm; 233c4762a1bSJed Brown DM dma, odm; 234c4762a1bSJed Brown Vec metric; 235c4762a1bSJed Brown size_t len; 236c4762a1bSJed Brown PetscErrorCode ierr; 237c4762a1bSJed Brown 238c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 239c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 240c4762a1bSJed Brown ierr = ProcessOptions(comm, &user);CHKERRQ(ierr); 241c4762a1bSJed Brown 242c4762a1bSJed Brown ierr = CreateMesh(comm, &user);CHKERRQ(ierr); 243c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) user.dm, "DMinit");CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = DMViewFromOptions(user.dm, NULL, "-init_dm_view");CHKERRQ(ierr); 245c4762a1bSJed Brown 246c4762a1bSJed Brown odm = user.dm; 247c4762a1bSJed Brown ierr = DMPlexDistributeOverlap(odm, 1, NULL, &user.dm);CHKERRQ(ierr); 248c4762a1bSJed Brown if (!user.dm) {user.dm = odm;} 249c4762a1bSJed Brown else {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 250c4762a1bSJed Brown ierr = ComputeMetric(user.dm, &user, &metric);CHKERRQ(ierr); 251c4762a1bSJed Brown ierr = PetscStrlen(user.bdLabel, &len);CHKERRQ(ierr); 252c4762a1bSJed Brown if (len) { 253c4762a1bSJed Brown ierr = DMCreateLabel(user.dm, user.bdLabel);CHKERRQ(ierr); 254c4762a1bSJed Brown ierr = DMGetLabel(user.dm, user.bdLabel, &bdLabel);CHKERRQ(ierr); 255c4762a1bSJed Brown } 256c4762a1bSJed Brown ierr = DMAdaptMetric(user.dm, metric, bdLabel, &dma);CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) dma, "DMadapt");CHKERRQ(ierr); 258c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject) dma, "adapt_");CHKERRQ(ierr); 259c4762a1bSJed Brown ierr = DMViewFromOptions(dma, NULL, "-dm_view");CHKERRQ(ierr); 260c4762a1bSJed Brown if (user.doL2) {ierr = TestL2Projection(user.dm, dma, &user);CHKERRQ(ierr);} 261c4762a1bSJed Brown ierr = DMDestroy(&dma);CHKERRQ(ierr); 262c4762a1bSJed Brown ierr = VecDestroy(&metric);CHKERRQ(ierr); 263c4762a1bSJed Brown ierr = DMDestroy(&user.dm);CHKERRQ(ierr); 264c4762a1bSJed Brown PetscFinalize(); 265c4762a1bSJed Brown return 0; 266c4762a1bSJed Brown } 267c4762a1bSJed Brown 268c4762a1bSJed Brown /*TEST 269c4762a1bSJed Brown 2708e3a2eefSMatthew G. Knepley build: 2718e3a2eefSMatthew G. Knepley requires: pragmatic 2728e3a2eefSMatthew G. Knepley 273c4762a1bSJed Brown test: 274c4762a1bSJed Brown suffix: 0 275c4762a1bSJed Brown args: -dim 2 -nbrVerEdge 5 -dm_plex_separate_marker 0 -met 2 -init_dm_view -adapt_dm_view 276c4762a1bSJed Brown test: 277c4762a1bSJed Brown suffix: 1 278c4762a1bSJed Brown args: -dim 2 -nbrVerEdge 5 -dm_plex_separate_marker 1 -bdLabel marker -met 2 -init_dm_view -adapt_dm_view 279c4762a1bSJed Brown test: 280c4762a1bSJed Brown suffix: 2 281c4762a1bSJed Brown args: -dim 3 -nbrVerEdge 5 -met 2 -init_dm_view -adapt_dm_view 282c4762a1bSJed Brown test: 283c4762a1bSJed Brown suffix: 3 284c4762a1bSJed Brown args: -dim 3 -nbrVerEdge 5 -bdLabel marker -met 2 -init_dm_view -adapt_dm_view 285c4762a1bSJed Brown test: 286c4762a1bSJed Brown suffix: 4 287c4762a1bSJed Brown nsize: 2 288c4762a1bSJed Brown args: -dim 2 -nbrVerEdge 3 -dm_plex_separate_marker 0 -met 2 -init_dm_view -adapt_dm_view 289c4762a1bSJed Brown test: 290c4762a1bSJed Brown suffix: 5 291c4762a1bSJed Brown nsize: 4 292c4762a1bSJed Brown args: -dim 2 -nbrVerEdge 3 -dm_plex_separate_marker 0 -met 2 -init_dm_view -adapt_dm_view 293c4762a1bSJed Brown test: 294c4762a1bSJed Brown suffix: 6 295c4762a1bSJed Brown nsize: 2 296*3e32e87aSJoe Wallwork args: -dim 3 -nbrVerEdge 10 -met 0 -hmin 0.01 -hmax 0.03 -init_dm_view -adapt_dm_view 297c4762a1bSJed Brown test: 298c4762a1bSJed Brown suffix: 7 299*3e32e87aSJoe Wallwork nsize: 2 300c4762a1bSJed Brown args: -dim 2 -nbrVerEdge 20 -dm_plex_separate_marker 0 -met 2 -hmax 0.5 -hmin 0.001 -init_dm_view -adapt_dm_view 301c4762a1bSJed Brown test: 302c4762a1bSJed Brown suffix: proj_0 303c4762a1bSJed Brown TODO: broken 304c4762a1bSJed Brown args: -dim 2 -nbrVerEdge 3 -dm_plex_separate_marker 0 -init_dm_view -adapt_dm_view -do_L2 -petscspace_degree 1 -petscfe_default_quadrature_order 1 -dm_plex_hash_location -pc_type lu 305c4762a1bSJed Brown test: 306c4762a1bSJed Brown suffix: proj_1 307c4762a1bSJed Brown TODO: broken 308c4762a1bSJed Brown args: -dim 2 -nbrVerEdge 5 -dm_plex_separate_marker 0 -init_dm_view -adapt_dm_view -do_L2 -petscspace_degree 2 -petscfe_default_quadrature_order 4 -dm_plex_hash_location -pc_type lu 309c4762a1bSJed Brown 310c4762a1bSJed Brown TEST*/ 311