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); 42*1e1ea65dSPierre 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 { 81c4762a1bSJed Brown DM cdm, mdm; 82c4762a1bSJed Brown PetscSection csec, msec; 83c4762a1bSJed Brown Vec coordinates; 84c4762a1bSJed Brown const PetscScalar *coords; 85c4762a1bSJed Brown PetscScalar *met; 86c4762a1bSJed Brown PetscReal h, *lambda, lbd, lmax; 87c4762a1bSJed Brown PetscInt pStart, pEnd, p, d; 88c4762a1bSJed Brown const PetscInt dim = user->dim, Nd = dim*dim; 89c4762a1bSJed Brown PetscErrorCode ierr; 90c4762a1bSJed Brown 91c4762a1bSJed Brown PetscFunctionBeginUser; 92c4762a1bSJed Brown ierr = PetscCalloc1(PetscMax(3, dim),&lambda);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 94c4762a1bSJed Brown ierr = DMClone(cdm, &mdm);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = DMGetLocalSection(cdm, &csec);CHKERRQ(ierr); 96c4762a1bSJed Brown 97c4762a1bSJed Brown ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &msec);CHKERRQ(ierr); 98c4762a1bSJed Brown ierr = PetscSectionSetNumFields(msec, 1);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = PetscSectionSetFieldComponents(msec, 0, Nd);CHKERRQ(ierr); 100c4762a1bSJed Brown ierr = PetscSectionGetChart(csec, &pStart, &pEnd);CHKERRQ(ierr); 101c4762a1bSJed Brown ierr = PetscSectionSetChart(msec, pStart, pEnd);CHKERRQ(ierr); 102c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 103c4762a1bSJed Brown ierr = PetscSectionSetDof(msec, p, Nd);CHKERRQ(ierr); 104c4762a1bSJed Brown ierr = PetscSectionSetFieldDof(msec, p, 0, Nd);CHKERRQ(ierr); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown ierr = PetscSectionSetUp(msec);CHKERRQ(ierr); 107c4762a1bSJed Brown ierr = DMSetLocalSection(mdm, msec);CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = PetscSectionDestroy(&msec);CHKERRQ(ierr); 109c4762a1bSJed Brown 110c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 111c4762a1bSJed Brown ierr = DMCreateLocalVector(mdm, metric);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = VecGetArray(*metric, &met);CHKERRQ(ierr); 114c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) { 115c4762a1bSJed Brown PetscScalar *pcoords; 116c4762a1bSJed Brown PetscScalar *pmet; 117c4762a1bSJed Brown 118c4762a1bSJed Brown ierr = DMPlexPointLocalRead(cdm, p, coords, &pcoords);CHKERRQ(ierr); 119c4762a1bSJed Brown switch (user->metOpt) { 120c4762a1bSJed Brown case 0: 121c4762a1bSJed Brown lbd = 1/(user->hmax*user->hmax); 122c4762a1bSJed Brown lambda[0] = lambda[1] = lambda[2] = lbd; 123c4762a1bSJed Brown break; 124c4762a1bSJed Brown case 1: 125c4762a1bSJed Brown h = user->hmax - (user->hmax-user->hmin)*PetscRealPart(pcoords[0]); 126c4762a1bSJed Brown h = h*h; 127c4762a1bSJed Brown lmax = 1/(user->hmax*user->hmax); 128c4762a1bSJed Brown lambda[0] = 1/h; 129c4762a1bSJed Brown lambda[1] = lmax; 130c4762a1bSJed Brown lambda[2] = lmax; 131c4762a1bSJed Brown break; 132c4762a1bSJed Brown case 2: 133c4762a1bSJed Brown h = user->hmax*PetscAbsReal(((PetscReal) 1.0)-PetscExpReal(-PetscAbsScalar(pcoords[0]-(PetscReal)0.5))) + user->hmin; 134c4762a1bSJed Brown lbd = 1/(h*h); 135c4762a1bSJed Brown lmax = 1/(user->hmax*user->hmax); 136c4762a1bSJed Brown lambda[0] = lbd; 137c4762a1bSJed Brown lambda[1] = lmax; 138c4762a1bSJed Brown lambda[2] = lmax; 139c4762a1bSJed Brown break; 140c4762a1bSJed Brown default: 141c4762a1bSJed Brown SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "metOpt = 0, 1 or 2, cannot be %d", user->metOpt); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown /* Only set the diagonal */ 144c4762a1bSJed Brown ierr = DMPlexPointLocalRef(mdm, p, met, &pmet);CHKERRQ(ierr); 145c4762a1bSJed Brown for (d = 0; d < dim; ++d) pmet[d*(dim+1)] = lambda[d]; 146c4762a1bSJed Brown } 147c4762a1bSJed Brown ierr = VecRestoreArray(*metric, &met);CHKERRQ(ierr); 148c4762a1bSJed Brown ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = DMDestroy(&mdm);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = PetscFree(lambda);CHKERRQ(ierr); 151c4762a1bSJed Brown PetscFunctionReturn(0); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 154c4762a1bSJed Brown static PetscErrorCode linear(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 155c4762a1bSJed Brown { 156c4762a1bSJed Brown u[0] = x[0] + x[1]; 157c4762a1bSJed Brown return 0; 158c4762a1bSJed Brown } 159c4762a1bSJed Brown 160c4762a1bSJed Brown static void identity(PetscInt dim, PetscInt Nf, PetscInt NfAux, 161c4762a1bSJed Brown const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], 162c4762a1bSJed Brown const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], 163c4762a1bSJed Brown PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[]) 164c4762a1bSJed Brown { 165c4762a1bSJed Brown g0[0] = 1.0; 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 168c4762a1bSJed Brown static PetscErrorCode TestL2Projection(DM dm, DM dma, AppCtx *user) 169c4762a1bSJed Brown { 170c4762a1bSJed Brown PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *); 171c4762a1bSJed Brown KSP ksp; 172c4762a1bSJed Brown PetscDS prob; 173c4762a1bSJed Brown PetscFE fe; 174c4762a1bSJed Brown Mat Interp, mass; 175c4762a1bSJed Brown Vec u, ua, scaling, ones, massLumped, rhs, uproj; 176c4762a1bSJed Brown PetscReal error; 177c4762a1bSJed Brown PetscInt dim; 178c4762a1bSJed Brown MPI_Comm comm; 179c4762a1bSJed Brown PetscErrorCode ierr; 180c4762a1bSJed Brown 181c4762a1bSJed Brown PetscFunctionBeginUser; 182c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 183c4762a1bSJed Brown ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 184c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, NULL, -1, &fe);CHKERRQ(ierr); 185c4762a1bSJed Brown ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = PetscDSSetDiscretization(prob, 0, (PetscObject) fe);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 188c4762a1bSJed Brown ierr = PetscFECreateDefault(comm, dim, 1, PETSC_TRUE, NULL, -1, &fe);CHKERRQ(ierr); 189c4762a1bSJed Brown ierr = DMGetDS(dma, &prob);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = PetscDSSetDiscretization(prob, 0, (PetscObject) fe);CHKERRQ(ierr); 191c4762a1bSJed Brown ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 192c4762a1bSJed Brown ierr = PetscDSSetJacobian(prob, 0, 0, identity, NULL, NULL, NULL);CHKERRQ(ierr); 193c4762a1bSJed Brown 194c4762a1bSJed Brown funcs[0] = linear; 195c4762a1bSJed Brown ierr = DMGetGlobalVector(dm, &u);CHKERRQ(ierr); 196c4762a1bSJed Brown ierr = DMGetGlobalVector(dma, &ua);CHKERRQ(ierr); 197c4762a1bSJed Brown ierr = DMGetGlobalVector(dma, &ones);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = DMGetGlobalVector(dma, &massLumped);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = DMGetGlobalVector(dma, &rhs);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = DMGetGlobalVector(dma, &uproj);CHKERRQ(ierr); 201c4762a1bSJed Brown ierr = DMProjectFunction(dm, 0.0, funcs, NULL, INSERT_VALUES, u);CHKERRQ(ierr); 202c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) u, "Original");CHKERRQ(ierr); 203c4762a1bSJed Brown ierr = VecViewFromOptions(u, NULL, "-orig_vec_view");CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = DMComputeL2Diff(dm, 0.0, funcs, NULL, u, &error);CHKERRQ(ierr); 205c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Original L2 Error: %g\n", (double) error);CHKERRQ(ierr); 206c4762a1bSJed Brown ierr = DMCreateInterpolation(dm, dma, &Interp, &scaling);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = MatInterpolate(Interp, u, ua);CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = MatDestroy(&Interp);CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = VecDestroy(&scaling);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) ua, "Interpolation");CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = VecViewFromOptions(ua, NULL, "-interp_vec_view");CHKERRQ(ierr); 212c4762a1bSJed Brown ierr = DMComputeL2Diff(dma, 0.0, funcs, NULL, ua, &error);CHKERRQ(ierr); 213c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Interpolated L2 Error: %g\n", (double) error);CHKERRQ(ierr); 214c4762a1bSJed Brown 215c4762a1bSJed Brown ierr = VecSet(ones, 1.0);CHKERRQ(ierr); 2168e3a2eefSMatthew G. Knepley ierr = DMSNESComputeJacobianAction(dma, ua, ones, massLumped, user);CHKERRQ(ierr); 217c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) massLumped, "Lumped mass");CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = VecViewFromOptions(massLumped, NULL, "-mass_vec_view");CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = DMCreateMassMatrix(dm, dma, &mass);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = MatMult(mass, u, rhs);CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = MatDestroy(&mass);CHKERRQ(ierr); 222c4762a1bSJed Brown ierr = VecViewFromOptions(rhs, NULL, "-lumped_rhs_view");CHKERRQ(ierr); 223c4762a1bSJed Brown ierr = VecPointwiseDivide(uproj, rhs, massLumped);CHKERRQ(ierr); 224c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) uproj, "Different Lumped Projection");CHKERRQ(ierr); 225c4762a1bSJed Brown ierr = VecViewFromOptions(uproj, NULL, "-lumped_rhs_vec_view");CHKERRQ(ierr); 226c4762a1bSJed Brown ierr = DMComputeL2Diff(dma, 0.0, funcs, NULL, uproj, &error);CHKERRQ(ierr); 227c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Lumped (rhs) L2 Error: %g\n", (double) error);CHKERRQ(ierr); 228c4762a1bSJed Brown 229c4762a1bSJed Brown ierr = DMCreateMatrix(dma, &mass);CHKERRQ(ierr); 230c4762a1bSJed Brown ierr = DMPlexSNESComputeJacobianFEM(dma, ua, mass, mass, user);CHKERRQ(ierr); 231c4762a1bSJed Brown ierr = MatViewFromOptions(mass, NULL, "-mass_mat_view");CHKERRQ(ierr); 232c4762a1bSJed Brown ierr = KSPCreate(PETSC_COMM_WORLD, &ksp);CHKERRQ(ierr); 233c4762a1bSJed Brown ierr = KSPSetOperators(ksp, mass, mass);CHKERRQ(ierr); 234c4762a1bSJed Brown ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = KSPSolve(ksp, rhs, uproj);CHKERRQ(ierr); 236c4762a1bSJed Brown ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 237c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) uproj, "Full Projection");CHKERRQ(ierr); 238c4762a1bSJed Brown ierr = VecViewFromOptions(uproj, NULL, "-proj_vec_view");CHKERRQ(ierr); 239c4762a1bSJed Brown ierr = DMComputeL2Diff(dma, 0.0, funcs, NULL, uproj, &error);CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "Projected L2 Error: %g\n", (double) error);CHKERRQ(ierr); 241c4762a1bSJed Brown ierr = MatDestroy(&mass);CHKERRQ(ierr); 242c4762a1bSJed Brown 243c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dm, &u);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dma, &ua);CHKERRQ(ierr); 245c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dma, &ones);CHKERRQ(ierr); 246c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dma, &massLumped);CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dma, &rhs);CHKERRQ(ierr); 248c4762a1bSJed Brown ierr = DMRestoreGlobalVector(dma, &uproj);CHKERRQ(ierr); 249c4762a1bSJed Brown PetscFunctionReturn(0); 250c4762a1bSJed Brown } 251c4762a1bSJed Brown 252c4762a1bSJed Brown int main (int argc, char * argv[]) { 253c4762a1bSJed Brown AppCtx user; /* user-defined work context */ 254c4762a1bSJed Brown DMLabel bdLabel = NULL; 255c4762a1bSJed Brown MPI_Comm comm; 256c4762a1bSJed Brown DM dma, odm; 257c4762a1bSJed Brown Vec metric; 258c4762a1bSJed Brown size_t len; 259c4762a1bSJed Brown PetscErrorCode ierr; 260c4762a1bSJed Brown 261c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 262c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 263c4762a1bSJed Brown ierr = ProcessOptions(comm, &user);CHKERRQ(ierr); 264c4762a1bSJed Brown 265c4762a1bSJed Brown ierr = CreateMesh(comm, &user);CHKERRQ(ierr); 266c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) user.dm, "DMinit");CHKERRQ(ierr); 267c4762a1bSJed Brown ierr = DMViewFromOptions(user.dm, NULL, "-init_dm_view");CHKERRQ(ierr); 268c4762a1bSJed Brown 269c4762a1bSJed Brown odm = user.dm; 270c4762a1bSJed Brown ierr = DMPlexDistributeOverlap(odm, 1, NULL, &user.dm);CHKERRQ(ierr); 271c4762a1bSJed Brown if (!user.dm) {user.dm = odm;} 272c4762a1bSJed Brown else {ierr = DMDestroy(&odm);CHKERRQ(ierr);} 273c4762a1bSJed Brown ierr = ComputeMetric(user.dm, &user, &metric);CHKERRQ(ierr); 274c4762a1bSJed Brown ierr = PetscStrlen(user.bdLabel, &len);CHKERRQ(ierr); 275c4762a1bSJed Brown if (len) { 276c4762a1bSJed Brown ierr = DMCreateLabel(user.dm, user.bdLabel);CHKERRQ(ierr); 277c4762a1bSJed Brown ierr = DMGetLabel(user.dm, user.bdLabel, &bdLabel);CHKERRQ(ierr); 278c4762a1bSJed Brown } 279c4762a1bSJed Brown ierr = DMAdaptMetric(user.dm, metric, bdLabel, &dma);CHKERRQ(ierr); 280c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject) dma, "DMadapt");CHKERRQ(ierr); 281c4762a1bSJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject) dma, "adapt_");CHKERRQ(ierr); 282c4762a1bSJed Brown ierr = DMViewFromOptions(dma, NULL, "-dm_view");CHKERRQ(ierr); 283c4762a1bSJed Brown if (user.doL2) {ierr = TestL2Projection(user.dm, dma, &user);CHKERRQ(ierr);} 284c4762a1bSJed Brown ierr = DMDestroy(&dma);CHKERRQ(ierr); 285c4762a1bSJed Brown ierr = VecDestroy(&metric);CHKERRQ(ierr); 286c4762a1bSJed Brown ierr = DMDestroy(&user.dm);CHKERRQ(ierr); 287c4762a1bSJed Brown PetscFinalize(); 288c4762a1bSJed Brown return 0; 289c4762a1bSJed Brown } 290c4762a1bSJed Brown 291c4762a1bSJed Brown /*TEST 292c4762a1bSJed Brown 2938e3a2eefSMatthew G. Knepley build: 2948e3a2eefSMatthew G. Knepley requires: pragmatic 2958e3a2eefSMatthew G. Knepley 296c4762a1bSJed Brown test: 297c4762a1bSJed Brown suffix: 0 298c4762a1bSJed Brown TODO: broken 299c4762a1bSJed Brown args: -dim 2 -nbrVerEdge 5 -dm_plex_separate_marker 0 -met 2 -init_dm_view -adapt_dm_view 300c4762a1bSJed Brown test: 301c4762a1bSJed Brown suffix: 1 302c4762a1bSJed Brown TODO: broken 303c4762a1bSJed Brown args: -dim 2 -nbrVerEdge 5 -dm_plex_separate_marker 1 -bdLabel marker -met 2 -init_dm_view -adapt_dm_view 304c4762a1bSJed Brown test: 305c4762a1bSJed Brown suffix: 2 306c4762a1bSJed Brown TODO: broken 307c4762a1bSJed Brown args: -dim 3 -nbrVerEdge 5 -met 2 -init_dm_view -adapt_dm_view 308c4762a1bSJed Brown test: 309c4762a1bSJed Brown suffix: 3 310c4762a1bSJed Brown TODO: broken 311c4762a1bSJed Brown args: -dim 3 -nbrVerEdge 5 -bdLabel marker -met 2 -init_dm_view -adapt_dm_view 312c4762a1bSJed Brown test: 313c4762a1bSJed Brown suffix: 4 314c4762a1bSJed Brown TODO: broken 315c4762a1bSJed Brown nsize: 2 316c4762a1bSJed Brown args: -dim 2 -nbrVerEdge 3 -dm_plex_separate_marker 0 -met 2 -init_dm_view -adapt_dm_view 317c4762a1bSJed Brown test: 318c4762a1bSJed Brown suffix: 5 319c4762a1bSJed Brown TODO: broken 320c4762a1bSJed Brown nsize: 4 321c4762a1bSJed Brown args: -dim 2 -nbrVerEdge 3 -dm_plex_separate_marker 0 -met 2 -init_dm_view -adapt_dm_view 322c4762a1bSJed Brown test: 323c4762a1bSJed Brown suffix: 6 324c4762a1bSJed Brown TODO: broken 325c4762a1bSJed Brown nsize: 2 326c4762a1bSJed Brown args: -dim 3 -nbrVerEdge 10 -dm_plex_separate_marker 0 -met 0 -hmin 0.01 -hmax 0.03 -init_dm_view -adapt_dm_view 327c4762a1bSJed Brown test: 328c4762a1bSJed Brown suffix: 7 329c4762a1bSJed Brown TODO: broken 330c4762a1bSJed Brown nsize: 5 331c4762a1bSJed 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 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown suffix: proj_0 334c4762a1bSJed Brown TODO: broken 335c4762a1bSJed 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 336c4762a1bSJed Brown test: 337c4762a1bSJed Brown suffix: proj_1 338c4762a1bSJed Brown TODO: broken 339c4762a1bSJed 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 340c4762a1bSJed Brown 341c4762a1bSJed Brown TEST*/ 342