1adcf2cb4SMatthew G. Knepley static char help[] = "TDycore Mesh Examples\n\n"; 2adcf2cb4SMatthew G. Knepley 3adcf2cb4SMatthew G. Knepley #include <petscdmplex.h> 4adcf2cb4SMatthew G. Knepley 5adcf2cb4SMatthew G. Knepley typedef struct { 6adcf2cb4SMatthew G. Knepley PetscBool adapt; /* Flag for adaptation of the surface mesh */ 7adcf2cb4SMatthew G. Knepley } AppCtx; 8adcf2cb4SMatthew G. Knepley 9adcf2cb4SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 10adcf2cb4SMatthew G. Knepley { 11adcf2cb4SMatthew G. Knepley PetscFunctionBeginUser; 12adcf2cb4SMatthew G. Knepley options->adapt = PETSC_FALSE; 13adcf2cb4SMatthew G. Knepley 14d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX"); 159566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-adapt", "Flag for adaptation of the surface mesh", "ex10.c", options->adapt, &options->adapt, NULL)); 16d0609cedSBarry Smith PetscOptionsEnd(); 17adcf2cb4SMatthew G. Knepley PetscFunctionReturn(0); 18adcf2cb4SMatthew G. Knepley } 19adcf2cb4SMatthew G. Knepley 20adcf2cb4SMatthew G. Knepley static PetscErrorCode CreateDomainLabel(DM dm) 21adcf2cb4SMatthew G. Knepley { 22adcf2cb4SMatthew G. Knepley DMLabel label; 23adcf2cb4SMatthew G. Knepley PetscInt cStart, cEnd, c; 24adcf2cb4SMatthew G. Knepley 25adcf2cb4SMatthew G. Knepley PetscFunctionBeginUser; 26*8fb5bd83SMatthew G. Knepley PetscCall(DMGetCoordinatesLocalSetUp(dm)); 279566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "Cell Sets")); 289566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Cell Sets", &label)); 299566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 30adcf2cb4SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 31adcf2cb4SMatthew G. Knepley PetscReal centroid[3], volume, x, y; 32adcf2cb4SMatthew G. Knepley 339566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &volume, centroid, NULL)); 34adcf2cb4SMatthew G. Knepley x = centroid[0]; y = centroid[1]; 35adcf2cb4SMatthew G. Knepley /* Headwaters are (0.0,0.25)--(0.1,0.75) */ 369566063dSJacob Faibussowitsch if ((x >= 0.0 && x < 0.1) && (y >= 0.25 && y <= 0.75)) {PetscCall(DMLabelSetValue(label, c, 1));continue;} 37adcf2cb4SMatthew G. Knepley /* River channel is (0.1,0.45)--(1.0,0.55) */ 389566063dSJacob Faibussowitsch if ((x >= 0.1 && x <= 1.0) && (y >= 0.45 && y <= 0.55)) {PetscCall(DMLabelSetValue(label, c, 2));continue;} 39adcf2cb4SMatthew G. Knepley } 40adcf2cb4SMatthew G. Knepley PetscFunctionReturn(0); 41adcf2cb4SMatthew G. Knepley } 42adcf2cb4SMatthew G. Knepley 43adcf2cb4SMatthew G. Knepley static PetscErrorCode AdaptMesh(DM *dm, AppCtx *ctx) 44adcf2cb4SMatthew G. Knepley { 45adcf2cb4SMatthew G. Knepley DM dmCur = *dm; 46adcf2cb4SMatthew G. Knepley DMLabel label; 47adcf2cb4SMatthew G. Knepley IS valueIS, vIS; 48adcf2cb4SMatthew G. Knepley PetscBool hasLabel; 49adcf2cb4SMatthew G. Knepley const PetscInt *values; 50adcf2cb4SMatthew G. Knepley PetscReal *volConst; /* Volume constraints for each label value */ 51adcf2cb4SMatthew G. Knepley PetscReal ratio; 52adcf2cb4SMatthew G. Knepley PetscInt dim, Nv, v, cStart, cEnd, c; 53adcf2cb4SMatthew G. Knepley PetscBool adapt = PETSC_TRUE; 54adcf2cb4SMatthew G. Knepley 55adcf2cb4SMatthew G. Knepley PetscFunctionBeginUser; 56adcf2cb4SMatthew G. Knepley if (!ctx->adapt) PetscFunctionReturn(0); 579566063dSJacob Faibussowitsch PetscCall(DMHasLabel(*dm, "Cell Sets", &hasLabel)); 589566063dSJacob Faibussowitsch if (!hasLabel) PetscCall(CreateDomainLabel(*dm)); 599566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 60adcf2cb4SMatthew G. Knepley ratio = PetscPowRealInt(0.5, dim); 61adcf2cb4SMatthew G. Knepley /* Get volume constraints */ 629566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "Cell Sets", &label)); 639566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(label, &vIS)); 649566063dSJacob Faibussowitsch PetscCall(ISDuplicate(vIS, &valueIS)); 659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vIS)); 66adcf2cb4SMatthew G. Knepley /* Sorting ruins the label */ 679566063dSJacob Faibussowitsch PetscCall(ISSort(valueIS)); 689566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(valueIS, &Nv)); 699566063dSJacob Faibussowitsch PetscCall(ISGetIndices(valueIS, &values)); 709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nv, &volConst)); 71adcf2cb4SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 72adcf2cb4SMatthew G. Knepley char opt[128]; 73adcf2cb4SMatthew G. Knepley 74adcf2cb4SMatthew G. Knepley volConst[v] = PETSC_MAX_REAL; 759566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(opt, 128, "-volume_constraint_%d", (int) values[v])); 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, opt, &volConst[v], NULL)); 77adcf2cb4SMatthew G. Knepley } 789566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(valueIS, &values)); 799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&valueIS)); 80adcf2cb4SMatthew G. Knepley /* Adapt mesh iteratively */ 81adcf2cb4SMatthew G. Knepley while (adapt) { 82adcf2cb4SMatthew G. Knepley DM dmAdapt; 83adcf2cb4SMatthew G. Knepley DMLabel adaptLabel; 84adcf2cb4SMatthew G. Knepley PetscInt nAdaptLoc[2], nAdapt[2]; 85adcf2cb4SMatthew G. Knepley 86adcf2cb4SMatthew G. Knepley adapt = PETSC_FALSE; 87adcf2cb4SMatthew G. Knepley nAdaptLoc[0] = nAdaptLoc[1] = 0; 88adcf2cb4SMatthew G. Knepley nAdapt[0] = nAdapt[1] = 0; 89adcf2cb4SMatthew G. Knepley /* Adaptation is not preserving the domain label */ 909566063dSJacob Faibussowitsch PetscCall(DMHasLabel(dmCur, "Cell Sets", &hasLabel)); 919566063dSJacob Faibussowitsch if (!hasLabel) PetscCall(CreateDomainLabel(dmCur)); 929566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dmCur, "Cell Sets", &label)); 939566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(label, &vIS)); 949566063dSJacob Faibussowitsch PetscCall(ISDuplicate(vIS, &valueIS)); 959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vIS)); 96adcf2cb4SMatthew G. Knepley /* Sorting directly the label's value IS would corrupt the label so we duplicate the IS first */ 979566063dSJacob Faibussowitsch PetscCall(ISSort(valueIS)); 989566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(valueIS, &Nv)); 999566063dSJacob Faibussowitsch PetscCall(ISGetIndices(valueIS, &values)); 100adcf2cb4SMatthew G. Knepley /* Construct adaptation label */ 1019566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 1029566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmCur, 0, &cStart, &cEnd)); 103adcf2cb4SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 104adcf2cb4SMatthew G. Knepley PetscReal volume, centroid[3]; 105adcf2cb4SMatthew G. Knepley PetscInt value, vidx; 106adcf2cb4SMatthew G. Knepley 1079566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dmCur, c, &volume, centroid, NULL)); 1089566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, c, &value)); 109adcf2cb4SMatthew G. Knepley if (value < 0) continue; 1109566063dSJacob Faibussowitsch PetscCall(PetscFindInt(value, Nv, values, &vidx)); 11163a3b9bcSJacob Faibussowitsch PetscCheck(vidx >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %" PetscInt_FMT " for cell %" PetscInt_FMT " does not exist in label", value, c); 1129566063dSJacob Faibussowitsch if (volume > volConst[vidx]) {PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE)); ++nAdaptLoc[0];} 1139566063dSJacob Faibussowitsch if (volume < volConst[vidx]*ratio) {PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_COARSEN)); ++nAdaptLoc[1];} 114adcf2cb4SMatthew G. Knepley } 1159566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(valueIS, &values)); 1169566063dSJacob Faibussowitsch PetscCall(ISDestroy(&valueIS)); 1179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&nAdaptLoc, &nAdapt, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) dmCur))); 118adcf2cb4SMatthew G. Knepley if (nAdapt[0]) { 11963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(dmCur, "Adapted mesh, marking %" PetscInt_FMT " cells for refinement, and %" PetscInt_FMT " cells for coarsening\n", nAdapt[0], nAdapt[1])); 1209566063dSJacob Faibussowitsch PetscCall(DMAdaptLabel(dmCur, adaptLabel, &dmAdapt)); 1219566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmCur)); 1229566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmAdapt, NULL, "-adapt_dm_view")); 123adcf2cb4SMatthew G. Knepley dmCur = dmAdapt; 124adcf2cb4SMatthew G. Knepley adapt = PETSC_TRUE; 125adcf2cb4SMatthew G. Knepley } 1269566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&adaptLabel)); 127adcf2cb4SMatthew G. Knepley } 1289566063dSJacob Faibussowitsch PetscCall(PetscFree(volConst)); 129adcf2cb4SMatthew G. Knepley *dm = dmCur; 130adcf2cb4SMatthew G. Knepley PetscFunctionReturn(0); 131adcf2cb4SMatthew G. Knepley } 132adcf2cb4SMatthew G. Knepley 133adcf2cb4SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 134adcf2cb4SMatthew G. Knepley { 135adcf2cb4SMatthew G. Knepley PetscInt dim; 136adcf2cb4SMatthew G. Knepley 137adcf2cb4SMatthew G. Knepley PetscFunctionBeginUser; 138adcf2cb4SMatthew G. Knepley /* Create top surface */ 1399566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1409566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1419566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, "init_")); 1429566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1439566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, NULL)); 144adcf2cb4SMatthew G. Knepley /* Adapt surface */ 1459566063dSJacob Faibussowitsch PetscCall(AdaptMesh(dm, user)); 146adcf2cb4SMatthew G. Knepley /* Extrude surface to get volume mesh */ 1479566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 1489566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinates(*dm)); 1499566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *dm, "Mesh")); 1509566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1519566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 152adcf2cb4SMatthew G. Knepley PetscFunctionReturn(0); 153adcf2cb4SMatthew G. Knepley } 154adcf2cb4SMatthew G. Knepley 155adcf2cb4SMatthew G. Knepley int main(int argc, char **argv) 156adcf2cb4SMatthew G. Knepley { 157adcf2cb4SMatthew G. Knepley DM dm; 158adcf2cb4SMatthew G. Knepley AppCtx user; 159adcf2cb4SMatthew G. Knepley 1609566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1619566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 1629566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 1639566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1649566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 165b122ec5aSJacob Faibussowitsch return 0; 166adcf2cb4SMatthew G. Knepley } 167adcf2cb4SMatthew G. Knepley 168adcf2cb4SMatthew G. Knepley /*TEST 169adcf2cb4SMatthew G. Knepley 170adcf2cb4SMatthew G. Knepley test: 171adcf2cb4SMatthew G. Knepley suffix: 0 172adcf2cb4SMatthew G. Knepley requires: triangle 173d410b0cfSMatthew G. Knepley args: -init_dm_plex_dim 2 -init_dm_plex_box_faces 1,1 -dm_extrude 1 -dm_view 174adcf2cb4SMatthew G. Knepley 175adcf2cb4SMatthew G. Knepley test: # Regularly refine the surface before extrusion 176adcf2cb4SMatthew G. Knepley suffix: 1 177adcf2cb4SMatthew G. Knepley requires: triangle 178d410b0cfSMatthew G. Knepley args: -init_dm_plex_dim 2 -init_dm_refine 2 -dm_extrude 1 -dm_view 179adcf2cb4SMatthew G. Knepley 180adcf2cb4SMatthew G. Knepley test: # Parallel run 181adcf2cb4SMatthew G. Knepley suffix: 2 182adcf2cb4SMatthew G. Knepley requires: triangle 183adcf2cb4SMatthew G. Knepley nsize: 5 184e600fa54SMatthew G. Knepley args: -init_dm_plex_dim 2 -init_dm_refine 3 -petscpartitioner_type simple -dm_extrude 3 -dm_view 185adcf2cb4SMatthew G. Knepley 186adcf2cb4SMatthew G. Knepley test: # adaptively refine the surface before extrusion 187adcf2cb4SMatthew G. Knepley suffix: 3 188adcf2cb4SMatthew G. Knepley requires: triangle 189d410b0cfSMatthew G. Knepley args: -init_dm_plex_dim 2 -init_dm_plex_box_faces 5,5 -adapt -volume_constraint_1 0.01 -volume_constraint_2 0.000625 -dm_extrude 10 190adcf2cb4SMatthew G. Knepley 191adcf2cb4SMatthew G. Knepley TEST*/ 192