1*adcf2cb4SMatthew G. Knepley static char help[] = "TDycore Mesh Examples\n\n"; 2*adcf2cb4SMatthew G. Knepley 3*adcf2cb4SMatthew G. Knepley #include <petscdmplex.h> 4*adcf2cb4SMatthew G. Knepley 5*adcf2cb4SMatthew G. Knepley typedef struct { 6*adcf2cb4SMatthew G. Knepley char filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */ 7*adcf2cb4SMatthew G. Knepley PetscBool adapt; /* Flag for adaptation of the surface mesh */ 8*adcf2cb4SMatthew G. Knepley PetscBool extrude; /* Flag for extrusion of the suraace mesh */ 9*adcf2cb4SMatthew G. Knepley } AppCtx; 10*adcf2cb4SMatthew G. Knepley 11*adcf2cb4SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 12*adcf2cb4SMatthew G. Knepley { 13*adcf2cb4SMatthew G. Knepley PetscErrorCode ierr; 14*adcf2cb4SMatthew G. Knepley 15*adcf2cb4SMatthew G. Knepley PetscFunctionBeginUser; 16*adcf2cb4SMatthew G. Knepley options->filename[0] = '\0'; 17*adcf2cb4SMatthew G. Knepley options->adapt = PETSC_FALSE; 18*adcf2cb4SMatthew G. Knepley options->extrude = PETSC_TRUE; 19*adcf2cb4SMatthew G. Knepley 20*adcf2cb4SMatthew G. Knepley ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr); 21*adcf2cb4SMatthew G. Knepley ierr = PetscOptionsString("-filename", "The mesh file", "ex10.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr); 22*adcf2cb4SMatthew G. Knepley ierr = PetscOptionsBool("-adapt", "Flag for adaptation of the surface mesh", "ex10.c", options->adapt, &options->adapt, NULL);CHKERRQ(ierr); 23*adcf2cb4SMatthew G. Knepley ierr = PetscOptionsBool("-extrude", "Flag for extrusion of the surface mesh", "ex10.c", options->extrude, &options->extrude, NULL);CHKERRQ(ierr); 24*adcf2cb4SMatthew G. Knepley ierr = PetscOptionsEnd(); 25*adcf2cb4SMatthew G. Knepley PetscFunctionReturn(0); 26*adcf2cb4SMatthew G. Knepley } 27*adcf2cb4SMatthew G. Knepley 28*adcf2cb4SMatthew G. Knepley static PetscErrorCode CreateDomainLabel(DM dm) 29*adcf2cb4SMatthew G. Knepley { 30*adcf2cb4SMatthew G. Knepley DMLabel label; 31*adcf2cb4SMatthew G. Knepley PetscInt cStart, cEnd, c; 32*adcf2cb4SMatthew G. Knepley PetscErrorCode ierr; 33*adcf2cb4SMatthew G. Knepley 34*adcf2cb4SMatthew G. Knepley PetscFunctionBeginUser; 35*adcf2cb4SMatthew G. Knepley ierr = DMCreateLabel(dm, "Cell Sets");CHKERRQ(ierr); 36*adcf2cb4SMatthew G. Knepley ierr = DMGetLabel(dm, "Cell Sets", &label);CHKERRQ(ierr); 37*adcf2cb4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 38*adcf2cb4SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 39*adcf2cb4SMatthew G. Knepley PetscReal centroid[3], volume, x, y; 40*adcf2cb4SMatthew G. Knepley 41*adcf2cb4SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dm, c, &volume, centroid, NULL);CHKERRQ(ierr); 42*adcf2cb4SMatthew G. Knepley x = centroid[0]; y = centroid[1]; 43*adcf2cb4SMatthew G. Knepley /* Headwaters are (0.0,0.25)--(0.1,0.75) */ 44*adcf2cb4SMatthew G. Knepley if ((x >= 0.0 && x < 0.1) && (y >= 0.25 && y <= 0.75)) {ierr = DMLabelSetValue(label, c, 1);CHKERRQ(ierr);continue;} 45*adcf2cb4SMatthew G. Knepley /* River channel is (0.1,0.45)--(1.0,0.55) */ 46*adcf2cb4SMatthew G. Knepley if ((x >= 0.1 && x <= 1.0) && (y >= 0.45 && y <= 0.55)) {ierr = DMLabelSetValue(label, c, 2);CHKERRQ(ierr);continue;} 47*adcf2cb4SMatthew G. Knepley } 48*adcf2cb4SMatthew G. Knepley PetscFunctionReturn(0); 49*adcf2cb4SMatthew G. Knepley } 50*adcf2cb4SMatthew G. Knepley 51*adcf2cb4SMatthew G. Knepley static PetscErrorCode AdaptMesh(DM *dm, AppCtx *ctx) 52*adcf2cb4SMatthew G. Knepley { 53*adcf2cb4SMatthew G. Knepley DM dmCur = *dm; 54*adcf2cb4SMatthew G. Knepley DMLabel label; 55*adcf2cb4SMatthew G. Knepley IS valueIS, vIS; 56*adcf2cb4SMatthew G. Knepley PetscBool hasLabel; 57*adcf2cb4SMatthew G. Knepley const PetscInt *values; 58*adcf2cb4SMatthew G. Knepley PetscReal *volConst; /* Volume constraints for each label value */ 59*adcf2cb4SMatthew G. Knepley PetscReal ratio; 60*adcf2cb4SMatthew G. Knepley PetscInt dim, Nv, v, cStart, cEnd, c; 61*adcf2cb4SMatthew G. Knepley PetscBool adapt = PETSC_TRUE; 62*adcf2cb4SMatthew G. Knepley PetscErrorCode ierr; 63*adcf2cb4SMatthew G. Knepley 64*adcf2cb4SMatthew G. Knepley PetscFunctionBeginUser; 65*adcf2cb4SMatthew G. Knepley if (!ctx->adapt) PetscFunctionReturn(0); 66*adcf2cb4SMatthew G. Knepley ierr = DMHasLabel(*dm, "Cell Sets", &hasLabel);CHKERRQ(ierr); 67*adcf2cb4SMatthew G. Knepley if (!hasLabel) {ierr = CreateDomainLabel(*dm);CHKERRQ(ierr);} 68*adcf2cb4SMatthew G. Knepley ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr); 69*adcf2cb4SMatthew G. Knepley ratio = PetscPowRealInt(0.5, dim); 70*adcf2cb4SMatthew G. Knepley /* Get volume constraints */ 71*adcf2cb4SMatthew G. Knepley ierr = DMGetLabel(*dm, "Cell Sets", &label);CHKERRQ(ierr); 72*adcf2cb4SMatthew G. Knepley ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr); 73*adcf2cb4SMatthew G. Knepley ierr = ISDuplicate(vIS, &valueIS);CHKERRQ(ierr); 74*adcf2cb4SMatthew G. Knepley ierr = ISDestroy(&vIS);CHKERRQ(ierr); 75*adcf2cb4SMatthew G. Knepley /* Sorting ruins the label */ 76*adcf2cb4SMatthew G. Knepley ierr = ISSort(valueIS);CHKERRQ(ierr); 77*adcf2cb4SMatthew G. Knepley ierr = ISGetLocalSize(valueIS, &Nv);CHKERRQ(ierr); 78*adcf2cb4SMatthew G. Knepley ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr); 79*adcf2cb4SMatthew G. Knepley ierr = PetscMalloc1(Nv, &volConst);CHKERRQ(ierr); 80*adcf2cb4SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 81*adcf2cb4SMatthew G. Knepley char opt[128]; 82*adcf2cb4SMatthew G. Knepley 83*adcf2cb4SMatthew G. Knepley volConst[v] = PETSC_MAX_REAL; 84*adcf2cb4SMatthew G. Knepley ierr = PetscSNPrintf(opt, 128, "-volume_constraint_%d", (int) values[v]);CHKERRQ(ierr); 85*adcf2cb4SMatthew G. Knepley ierr = PetscOptionsGetReal(NULL, NULL, opt, &volConst[v], NULL);CHKERRQ(ierr); 86*adcf2cb4SMatthew G. Knepley } 87*adcf2cb4SMatthew G. Knepley ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr); 88*adcf2cb4SMatthew G. Knepley ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 89*adcf2cb4SMatthew G. Knepley /* Adapt mesh iteratively */ 90*adcf2cb4SMatthew G. Knepley while (adapt) { 91*adcf2cb4SMatthew G. Knepley DM dmAdapt; 92*adcf2cb4SMatthew G. Knepley DMLabel adaptLabel; 93*adcf2cb4SMatthew G. Knepley PetscInt nAdaptLoc[2], nAdapt[2]; 94*adcf2cb4SMatthew G. Knepley 95*adcf2cb4SMatthew G. Knepley adapt = PETSC_FALSE; 96*adcf2cb4SMatthew G. Knepley nAdaptLoc[0] = nAdaptLoc[1] = 0; 97*adcf2cb4SMatthew G. Knepley nAdapt[0] = nAdapt[1] = 0; 98*adcf2cb4SMatthew G. Knepley /* Adaptation is not preserving the domain label */ 99*adcf2cb4SMatthew G. Knepley ierr = DMHasLabel(dmCur, "Cell Sets", &hasLabel);CHKERRQ(ierr); 100*adcf2cb4SMatthew G. Knepley if (!hasLabel) {ierr = CreateDomainLabel(dmCur);CHKERRQ(ierr);} 101*adcf2cb4SMatthew G. Knepley ierr = DMGetLabel(dmCur, "Cell Sets", &label);CHKERRQ(ierr); 102*adcf2cb4SMatthew G. Knepley ierr = DMLabelGetValueIS(label, &vIS);CHKERRQ(ierr); 103*adcf2cb4SMatthew G. Knepley ierr = ISDuplicate(vIS, &valueIS);CHKERRQ(ierr); 104*adcf2cb4SMatthew G. Knepley ierr = ISDestroy(&vIS);CHKERRQ(ierr); 105*adcf2cb4SMatthew G. Knepley /* Sorting directly the label's value IS would corrupt the label so we duplicate the IS first */ 106*adcf2cb4SMatthew G. Knepley ierr = ISSort(valueIS);CHKERRQ(ierr); 107*adcf2cb4SMatthew G. Knepley ierr = ISGetLocalSize(valueIS, &Nv);CHKERRQ(ierr); 108*adcf2cb4SMatthew G. Knepley ierr = ISGetIndices(valueIS, &values);CHKERRQ(ierr); 109*adcf2cb4SMatthew G. Knepley /* Construct adaptation label */ 110*adcf2cb4SMatthew G. Knepley ierr = DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel);CHKERRQ(ierr); 111*adcf2cb4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dmCur, 0, &cStart, &cEnd);CHKERRQ(ierr); 112*adcf2cb4SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 113*adcf2cb4SMatthew G. Knepley PetscReal volume, centroid[3]; 114*adcf2cb4SMatthew G. Knepley PetscInt value, vidx; 115*adcf2cb4SMatthew G. Knepley 116*adcf2cb4SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dmCur, c, &volume, centroid, NULL);CHKERRQ(ierr); 117*adcf2cb4SMatthew G. Knepley ierr = DMLabelGetValue(label, c, &value);CHKERRQ(ierr); 118*adcf2cb4SMatthew G. Knepley if (value < 0) continue; 119*adcf2cb4SMatthew G. Knepley ierr = PetscFindInt(value, Nv, values, &vidx);CHKERRQ(ierr); 120*adcf2cb4SMatthew G. Knepley if (vidx < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %D for cell %D does not exist in label", value, c); 121*adcf2cb4SMatthew G. Knepley if (volume > volConst[vidx]) {ierr = DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE);CHKERRQ(ierr); ++nAdaptLoc[0];} 122*adcf2cb4SMatthew G. Knepley if (volume < volConst[vidx]*ratio) {ierr = DMLabelSetValue(adaptLabel, c, DM_ADAPT_COARSEN);CHKERRQ(ierr); ++nAdaptLoc[1];} 123*adcf2cb4SMatthew G. Knepley } 124*adcf2cb4SMatthew G. Knepley ierr = ISRestoreIndices(valueIS, &values);CHKERRQ(ierr); 125*adcf2cb4SMatthew G. Knepley ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 126*adcf2cb4SMatthew G. Knepley ierr = MPI_Allreduce(&nAdaptLoc, &nAdapt, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) dmCur));CHKERRQ(ierr); 127*adcf2cb4SMatthew G. Knepley if (nAdapt[0]) { 128*adcf2cb4SMatthew G. Knepley ierr = PetscInfo2(dmCur, "Adapted mesh, marking %D cells for refinement, and %D cells for coarsening\n", nAdapt[0], nAdapt[1]);CHKERRQ(ierr); 129*adcf2cb4SMatthew G. Knepley ierr = DMAdaptLabel(dmCur, adaptLabel, &dmAdapt);CHKERRQ(ierr); 130*adcf2cb4SMatthew G. Knepley ierr = DMDestroy(&dmCur); 131*adcf2cb4SMatthew G. Knepley ierr = DMViewFromOptions(dmAdapt, NULL, "-adapt_dm_view");CHKERRQ(ierr); 132*adcf2cb4SMatthew G. Knepley dmCur = dmAdapt; 133*adcf2cb4SMatthew G. Knepley adapt = PETSC_TRUE; 134*adcf2cb4SMatthew G. Knepley } 135*adcf2cb4SMatthew G. Knepley ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr); 136*adcf2cb4SMatthew G. Knepley } 137*adcf2cb4SMatthew G. Knepley ierr = PetscFree(volConst);CHKERRQ(ierr); 138*adcf2cb4SMatthew G. Knepley *dm = dmCur; 139*adcf2cb4SMatthew G. Knepley PetscFunctionReturn(0); 140*adcf2cb4SMatthew G. Knepley } 141*adcf2cb4SMatthew G. Knepley 142*adcf2cb4SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 143*adcf2cb4SMatthew G. Knepley { 144*adcf2cb4SMatthew G. Knepley PetscInt dim; 145*adcf2cb4SMatthew G. Knepley size_t len; 146*adcf2cb4SMatthew G. Knepley PetscErrorCode ierr; 147*adcf2cb4SMatthew G. Knepley 148*adcf2cb4SMatthew G. Knepley PetscFunctionBeginUser; 149*adcf2cb4SMatthew G. Knepley /* Create top surface */ 150*adcf2cb4SMatthew G. Knepley ierr = PetscStrlen(user->filename, &len);CHKERRQ(ierr); 151*adcf2cb4SMatthew G. Knepley if (len) {ierr = DMPlexCreateFromFile(comm, user->filename, PETSC_TRUE, dm);CHKERRQ(ierr);} 152*adcf2cb4SMatthew G. Knepley else {ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);} 153*adcf2cb4SMatthew G. Knepley /* Adapt surface */ 154*adcf2cb4SMatthew G. Knepley ierr = AdaptMesh(dm, user);CHKERRQ(ierr); 155*adcf2cb4SMatthew G. Knepley /* Extrude surface to get volume mesh */ 156*adcf2cb4SMatthew G. Knepley ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr); 157*adcf2cb4SMatthew G. Knepley if (dim < 3 && user->extrude) { 158*adcf2cb4SMatthew G. Knepley DM edm; 159*adcf2cb4SMatthew G. Knepley 160*adcf2cb4SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, "srf_");CHKERRQ(ierr); 161*adcf2cb4SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 162*adcf2cb4SMatthew G. Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 163*adcf2cb4SMatthew G. Knepley ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, NULL);CHKERRQ(ierr); 164*adcf2cb4SMatthew G. Knepley ierr = DMPlexExtrude(*dm, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_TRUE, NULL, PETSC_TRUE, &edm);CHKERRQ(ierr); 165*adcf2cb4SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 166*adcf2cb4SMatthew G. Knepley *dm = edm; 167*adcf2cb4SMatthew G. Knepley } 168*adcf2cb4SMatthew G. Knepley ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr); 169*adcf2cb4SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr); 170*adcf2cb4SMatthew G. Knepley ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 171*adcf2cb4SMatthew G. Knepley ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 172*adcf2cb4SMatthew G. Knepley PetscFunctionReturn(0); 173*adcf2cb4SMatthew G. Knepley } 174*adcf2cb4SMatthew G. Knepley 175*adcf2cb4SMatthew G. Knepley int main(int argc, char **argv) 176*adcf2cb4SMatthew G. Knepley { 177*adcf2cb4SMatthew G. Knepley DM dm; 178*adcf2cb4SMatthew G. Knepley AppCtx user; 179*adcf2cb4SMatthew G. Knepley PetscErrorCode ierr; 180*adcf2cb4SMatthew G. Knepley 181*adcf2cb4SMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 182*adcf2cb4SMatthew G. Knepley ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr); 183*adcf2cb4SMatthew G. Knepley ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr); 184*adcf2cb4SMatthew G. Knepley ierr = DMDestroy(&dm);CHKERRQ(ierr); 185*adcf2cb4SMatthew G. Knepley ierr = PetscFinalize(); 186*adcf2cb4SMatthew G. Knepley return ierr; 187*adcf2cb4SMatthew G. Knepley } 188*adcf2cb4SMatthew G. Knepley 189*adcf2cb4SMatthew G. Knepley /*TEST 190*adcf2cb4SMatthew G. Knepley 191*adcf2cb4SMatthew G. Knepley test: 192*adcf2cb4SMatthew G. Knepley suffix: 0 193*adcf2cb4SMatthew G. Knepley requires: triangle 194*adcf2cb4SMatthew G. Knepley args: -dm_plex_box_dim 2 -dm_plex_box_faces 1,1 -dm_view 195*adcf2cb4SMatthew G. Knepley 196*adcf2cb4SMatthew G. Knepley test: # Regularly refine the surface before extrusion 197*adcf2cb4SMatthew G. Knepley suffix: 1 198*adcf2cb4SMatthew G. Knepley requires: triangle 199*adcf2cb4SMatthew G. Knepley args: -dm_plex_box_dim 2 -srf_dm_refine 2 -dm_view 200*adcf2cb4SMatthew G. Knepley 201*adcf2cb4SMatthew G. Knepley test: # Parallel run 202*adcf2cb4SMatthew G. Knepley suffix: 2 203*adcf2cb4SMatthew G. Knepley requires: triangle 204*adcf2cb4SMatthew G. Knepley nsize: 5 205*adcf2cb4SMatthew G. Knepley args: -dm_plex_box_dim 2 -srf_dm_refine 3 -petscpartitioner_type simple -dm_distribute -dm_plex_extrude_layers 3 -dm_view 206*adcf2cb4SMatthew G. Knepley 207*adcf2cb4SMatthew G. Knepley test: # adaptively refine the surface before extrusion 208*adcf2cb4SMatthew G. Knepley suffix: 3 209*adcf2cb4SMatthew G. Knepley requires: triangle 210*adcf2cb4SMatthew G. Knepley args: -dm_plex_box_dim 2 -dm_plex_box_faces 5,5 -adapt -volume_constraint_1 0.01 -volume_constraint_2 0.000625 -dm_plex_extrude_layers 10 211*adcf2cb4SMatthew G. Knepley 212*adcf2cb4SMatthew G. Knepley TEST*/ 213