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 14*d0609cedSBarry 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)); 16*d0609cedSBarry 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; 269566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "Cell Sets")); 279566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Cell Sets", &label)); 289566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 29adcf2cb4SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 30adcf2cb4SMatthew G. Knepley PetscReal centroid[3], volume, x, y; 31adcf2cb4SMatthew G. Knepley 329566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dm, c, &volume, centroid, NULL)); 33adcf2cb4SMatthew G. Knepley x = centroid[0]; y = centroid[1]; 34adcf2cb4SMatthew G. Knepley /* Headwaters are (0.0,0.25)--(0.1,0.75) */ 359566063dSJacob Faibussowitsch if ((x >= 0.0 && x < 0.1) && (y >= 0.25 && y <= 0.75)) {PetscCall(DMLabelSetValue(label, c, 1));continue;} 36adcf2cb4SMatthew G. Knepley /* River channel is (0.1,0.45)--(1.0,0.55) */ 379566063dSJacob Faibussowitsch if ((x >= 0.1 && x <= 1.0) && (y >= 0.45 && y <= 0.55)) {PetscCall(DMLabelSetValue(label, c, 2));continue;} 38adcf2cb4SMatthew G. Knepley } 39adcf2cb4SMatthew G. Knepley PetscFunctionReturn(0); 40adcf2cb4SMatthew G. Knepley } 41adcf2cb4SMatthew G. Knepley 42adcf2cb4SMatthew G. Knepley static PetscErrorCode AdaptMesh(DM *dm, AppCtx *ctx) 43adcf2cb4SMatthew G. Knepley { 44adcf2cb4SMatthew G. Knepley DM dmCur = *dm; 45adcf2cb4SMatthew G. Knepley DMLabel label; 46adcf2cb4SMatthew G. Knepley IS valueIS, vIS; 47adcf2cb4SMatthew G. Knepley PetscBool hasLabel; 48adcf2cb4SMatthew G. Knepley const PetscInt *values; 49adcf2cb4SMatthew G. Knepley PetscReal *volConst; /* Volume constraints for each label value */ 50adcf2cb4SMatthew G. Knepley PetscReal ratio; 51adcf2cb4SMatthew G. Knepley PetscInt dim, Nv, v, cStart, cEnd, c; 52adcf2cb4SMatthew G. Knepley PetscBool adapt = PETSC_TRUE; 53adcf2cb4SMatthew G. Knepley 54adcf2cb4SMatthew G. Knepley PetscFunctionBeginUser; 55adcf2cb4SMatthew G. Knepley if (!ctx->adapt) PetscFunctionReturn(0); 569566063dSJacob Faibussowitsch PetscCall(DMHasLabel(*dm, "Cell Sets", &hasLabel)); 579566063dSJacob Faibussowitsch if (!hasLabel) PetscCall(CreateDomainLabel(*dm)); 589566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 59adcf2cb4SMatthew G. Knepley ratio = PetscPowRealInt(0.5, dim); 60adcf2cb4SMatthew G. Knepley /* Get volume constraints */ 619566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "Cell Sets", &label)); 629566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(label, &vIS)); 639566063dSJacob Faibussowitsch PetscCall(ISDuplicate(vIS, &valueIS)); 649566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vIS)); 65adcf2cb4SMatthew G. Knepley /* Sorting ruins the label */ 669566063dSJacob Faibussowitsch PetscCall(ISSort(valueIS)); 679566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(valueIS, &Nv)); 689566063dSJacob Faibussowitsch PetscCall(ISGetIndices(valueIS, &values)); 699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nv, &volConst)); 70adcf2cb4SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 71adcf2cb4SMatthew G. Knepley char opt[128]; 72adcf2cb4SMatthew G. Knepley 73adcf2cb4SMatthew G. Knepley volConst[v] = PETSC_MAX_REAL; 749566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(opt, 128, "-volume_constraint_%d", (int) values[v])); 759566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, opt, &volConst[v], NULL)); 76adcf2cb4SMatthew G. Knepley } 779566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(valueIS, &values)); 789566063dSJacob Faibussowitsch PetscCall(ISDestroy(&valueIS)); 79adcf2cb4SMatthew G. Knepley /* Adapt mesh iteratively */ 80adcf2cb4SMatthew G. Knepley while (adapt) { 81adcf2cb4SMatthew G. Knepley DM dmAdapt; 82adcf2cb4SMatthew G. Knepley DMLabel adaptLabel; 83adcf2cb4SMatthew G. Knepley PetscInt nAdaptLoc[2], nAdapt[2]; 84adcf2cb4SMatthew G. Knepley 85adcf2cb4SMatthew G. Knepley adapt = PETSC_FALSE; 86adcf2cb4SMatthew G. Knepley nAdaptLoc[0] = nAdaptLoc[1] = 0; 87adcf2cb4SMatthew G. Knepley nAdapt[0] = nAdapt[1] = 0; 88adcf2cb4SMatthew G. Knepley /* Adaptation is not preserving the domain label */ 899566063dSJacob Faibussowitsch PetscCall(DMHasLabel(dmCur, "Cell Sets", &hasLabel)); 909566063dSJacob Faibussowitsch if (!hasLabel) PetscCall(CreateDomainLabel(dmCur)); 919566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dmCur, "Cell Sets", &label)); 929566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(label, &vIS)); 939566063dSJacob Faibussowitsch PetscCall(ISDuplicate(vIS, &valueIS)); 949566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vIS)); 95adcf2cb4SMatthew G. Knepley /* Sorting directly the label's value IS would corrupt the label so we duplicate the IS first */ 969566063dSJacob Faibussowitsch PetscCall(ISSort(valueIS)); 979566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(valueIS, &Nv)); 989566063dSJacob Faibussowitsch PetscCall(ISGetIndices(valueIS, &values)); 99adcf2cb4SMatthew G. Knepley /* Construct adaptation label */ 1009566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 1019566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmCur, 0, &cStart, &cEnd)); 102adcf2cb4SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 103adcf2cb4SMatthew G. Knepley PetscReal volume, centroid[3]; 104adcf2cb4SMatthew G. Knepley PetscInt value, vidx; 105adcf2cb4SMatthew G. Knepley 1069566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dmCur, c, &volume, centroid, NULL)); 1079566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, c, &value)); 108adcf2cb4SMatthew G. Knepley if (value < 0) continue; 1099566063dSJacob Faibussowitsch PetscCall(PetscFindInt(value, Nv, values, &vidx)); 11008401ef6SPierre Jolivet PetscCheck(vidx >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Value %D for cell %D does not exist in label", value, c); 1119566063dSJacob Faibussowitsch if (volume > volConst[vidx]) {PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE)); ++nAdaptLoc[0];} 1129566063dSJacob Faibussowitsch if (volume < volConst[vidx]*ratio) {PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_COARSEN)); ++nAdaptLoc[1];} 113adcf2cb4SMatthew G. Knepley } 1149566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(valueIS, &values)); 1159566063dSJacob Faibussowitsch PetscCall(ISDestroy(&valueIS)); 1169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&nAdaptLoc, &nAdapt, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) dmCur))); 117adcf2cb4SMatthew G. Knepley if (nAdapt[0]) { 1189566063dSJacob Faibussowitsch PetscCall(PetscInfo(dmCur, "Adapted mesh, marking %D cells for refinement, and %D cells for coarsening\n", nAdapt[0], nAdapt[1])); 1199566063dSJacob Faibussowitsch PetscCall(DMAdaptLabel(dmCur, adaptLabel, &dmAdapt)); 1209566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmCur)); 1219566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmAdapt, NULL, "-adapt_dm_view")); 122adcf2cb4SMatthew G. Knepley dmCur = dmAdapt; 123adcf2cb4SMatthew G. Knepley adapt = PETSC_TRUE; 124adcf2cb4SMatthew G. Knepley } 1259566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&adaptLabel)); 126adcf2cb4SMatthew G. Knepley } 1279566063dSJacob Faibussowitsch PetscCall(PetscFree(volConst)); 128adcf2cb4SMatthew G. Knepley *dm = dmCur; 129adcf2cb4SMatthew G. Knepley PetscFunctionReturn(0); 130adcf2cb4SMatthew G. Knepley } 131adcf2cb4SMatthew G. Knepley 132adcf2cb4SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 133adcf2cb4SMatthew G. Knepley { 134adcf2cb4SMatthew G. Knepley PetscInt dim; 135adcf2cb4SMatthew G. Knepley 136adcf2cb4SMatthew G. Knepley PetscFunctionBeginUser; 137adcf2cb4SMatthew G. Knepley /* Create top surface */ 1389566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1399566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1409566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, "init_")); 1419566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1429566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, NULL)); 143adcf2cb4SMatthew G. Knepley /* Adapt surface */ 1449566063dSJacob Faibussowitsch PetscCall(AdaptMesh(dm, user)); 145adcf2cb4SMatthew G. Knepley /* Extrude surface to get volume mesh */ 1469566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 1479566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinates(*dm)); 1489566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *dm, "Mesh")); 1499566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1509566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 151adcf2cb4SMatthew G. Knepley PetscFunctionReturn(0); 152adcf2cb4SMatthew G. Knepley } 153adcf2cb4SMatthew G. Knepley 154adcf2cb4SMatthew G. Knepley int main(int argc, char **argv) 155adcf2cb4SMatthew G. Knepley { 156adcf2cb4SMatthew G. Knepley DM dm; 157adcf2cb4SMatthew G. Knepley AppCtx user; 158adcf2cb4SMatthew G. Knepley 1599566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1609566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 1619566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 1629566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1639566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 164b122ec5aSJacob Faibussowitsch return 0; 165adcf2cb4SMatthew G. Knepley } 166adcf2cb4SMatthew G. Knepley 167adcf2cb4SMatthew G. Knepley /*TEST 168adcf2cb4SMatthew G. Knepley 169adcf2cb4SMatthew G. Knepley test: 170adcf2cb4SMatthew G. Knepley suffix: 0 171adcf2cb4SMatthew G. Knepley requires: triangle 172d410b0cfSMatthew G. Knepley args: -init_dm_plex_dim 2 -init_dm_plex_box_faces 1,1 -dm_extrude 1 -dm_view 173adcf2cb4SMatthew G. Knepley 174adcf2cb4SMatthew G. Knepley test: # Regularly refine the surface before extrusion 175adcf2cb4SMatthew G. Knepley suffix: 1 176adcf2cb4SMatthew G. Knepley requires: triangle 177d410b0cfSMatthew G. Knepley args: -init_dm_plex_dim 2 -init_dm_refine 2 -dm_extrude 1 -dm_view 178adcf2cb4SMatthew G. Knepley 179adcf2cb4SMatthew G. Knepley test: # Parallel run 180adcf2cb4SMatthew G. Knepley suffix: 2 181adcf2cb4SMatthew G. Knepley requires: triangle 182adcf2cb4SMatthew G. Knepley nsize: 5 183e600fa54SMatthew G. Knepley args: -init_dm_plex_dim 2 -init_dm_refine 3 -petscpartitioner_type simple -dm_extrude 3 -dm_view 184adcf2cb4SMatthew G. Knepley 185adcf2cb4SMatthew G. Knepley test: # adaptively refine the surface before extrusion 186adcf2cb4SMatthew G. Knepley suffix: 3 187adcf2cb4SMatthew G. Knepley requires: triangle 188d410b0cfSMatthew 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 189adcf2cb4SMatthew G. Knepley 190adcf2cb4SMatthew G. Knepley TEST*/ 191