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 9d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 10d71ae5a4SJacob Faibussowitsch { 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(); 173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18adcf2cb4SMatthew G. Knepley } 19adcf2cb4SMatthew G. Knepley 20d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateDomainLabel(DM dm) 21d71ae5a4SJacob Faibussowitsch { 22adcf2cb4SMatthew G. Knepley DMLabel label; 23adcf2cb4SMatthew G. Knepley PetscInt cStart, cEnd, c; 24adcf2cb4SMatthew G. Knepley 25adcf2cb4SMatthew G. Knepley PetscFunctionBeginUser; 268fb5bd83SMatthew 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)); 349371c9d4SSatish Balay x = centroid[0]; 359371c9d4SSatish Balay y = centroid[1]; 36adcf2cb4SMatthew G. Knepley /* Headwaters are (0.0,0.25)--(0.1,0.75) */ 379371c9d4SSatish Balay if ((x >= 0.0 && x < 0.1) && (y >= 0.25 && y <= 0.75)) { 389371c9d4SSatish Balay PetscCall(DMLabelSetValue(label, c, 1)); 399371c9d4SSatish Balay continue; 409371c9d4SSatish Balay } 41adcf2cb4SMatthew G. Knepley /* River channel is (0.1,0.45)--(1.0,0.55) */ 429371c9d4SSatish Balay if ((x >= 0.1 && x <= 1.0) && (y >= 0.45 && y <= 0.55)) { 439371c9d4SSatish Balay PetscCall(DMLabelSetValue(label, c, 2)); 449371c9d4SSatish Balay continue; 459371c9d4SSatish Balay } 46adcf2cb4SMatthew G. Knepley } 473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48adcf2cb4SMatthew G. Knepley } 49adcf2cb4SMatthew G. Knepley 50d71ae5a4SJacob Faibussowitsch static PetscErrorCode AdaptMesh(DM *dm, AppCtx *ctx) 51d71ae5a4SJacob Faibussowitsch { 52adcf2cb4SMatthew G. Knepley DM dmCur = *dm; 53adcf2cb4SMatthew G. Knepley DMLabel label; 54adcf2cb4SMatthew G. Knepley IS valueIS, vIS; 55adcf2cb4SMatthew G. Knepley PetscBool hasLabel; 56adcf2cb4SMatthew G. Knepley const PetscInt *values; 57adcf2cb4SMatthew G. Knepley PetscReal *volConst; /* Volume constraints for each label value */ 58adcf2cb4SMatthew G. Knepley PetscReal ratio; 59adcf2cb4SMatthew G. Knepley PetscInt dim, Nv, v, cStart, cEnd, c; 60adcf2cb4SMatthew G. Knepley PetscBool adapt = PETSC_TRUE; 61adcf2cb4SMatthew G. Knepley 62adcf2cb4SMatthew G. Knepley PetscFunctionBeginUser; 633ba16761SJacob Faibussowitsch if (!ctx->adapt) PetscFunctionReturn(PETSC_SUCCESS); 649566063dSJacob Faibussowitsch PetscCall(DMHasLabel(*dm, "Cell Sets", &hasLabel)); 659566063dSJacob Faibussowitsch if (!hasLabel) PetscCall(CreateDomainLabel(*dm)); 669566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 67adcf2cb4SMatthew G. Knepley ratio = PetscPowRealInt(0.5, dim); 68adcf2cb4SMatthew G. Knepley /* Get volume constraints */ 699566063dSJacob Faibussowitsch PetscCall(DMGetLabel(*dm, "Cell Sets", &label)); 709566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(label, &vIS)); 719566063dSJacob Faibussowitsch PetscCall(ISDuplicate(vIS, &valueIS)); 729566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vIS)); 73adcf2cb4SMatthew G. Knepley /* Sorting ruins the label */ 749566063dSJacob Faibussowitsch PetscCall(ISSort(valueIS)); 759566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(valueIS, &Nv)); 769566063dSJacob Faibussowitsch PetscCall(ISGetIndices(valueIS, &values)); 779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nv, &volConst)); 78adcf2cb4SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 79adcf2cb4SMatthew G. Knepley char opt[128]; 80adcf2cb4SMatthew G. Knepley 81adcf2cb4SMatthew G. Knepley volConst[v] = PETSC_MAX_REAL; 829566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(opt, 128, "-volume_constraint_%d", (int)values[v])); 839566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, opt, &volConst[v], NULL)); 84adcf2cb4SMatthew G. Knepley } 859566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(valueIS, &values)); 869566063dSJacob Faibussowitsch PetscCall(ISDestroy(&valueIS)); 87adcf2cb4SMatthew G. Knepley /* Adapt mesh iteratively */ 88adcf2cb4SMatthew G. Knepley while (adapt) { 89adcf2cb4SMatthew G. Knepley DM dmAdapt; 90adcf2cb4SMatthew G. Knepley DMLabel adaptLabel; 91adcf2cb4SMatthew G. Knepley PetscInt nAdaptLoc[2], nAdapt[2]; 92adcf2cb4SMatthew G. Knepley 93adcf2cb4SMatthew G. Knepley adapt = PETSC_FALSE; 94adcf2cb4SMatthew G. Knepley nAdaptLoc[0] = nAdaptLoc[1] = 0; 95adcf2cb4SMatthew G. Knepley nAdapt[0] = nAdapt[1] = 0; 96adcf2cb4SMatthew G. Knepley /* Adaptation is not preserving the domain label */ 979566063dSJacob Faibussowitsch PetscCall(DMHasLabel(dmCur, "Cell Sets", &hasLabel)); 989566063dSJacob Faibussowitsch if (!hasLabel) PetscCall(CreateDomainLabel(dmCur)); 999566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dmCur, "Cell Sets", &label)); 1009566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(label, &vIS)); 1019566063dSJacob Faibussowitsch PetscCall(ISDuplicate(vIS, &valueIS)); 1029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vIS)); 103adcf2cb4SMatthew G. Knepley /* Sorting directly the label's value IS would corrupt the label so we duplicate the IS first */ 1049566063dSJacob Faibussowitsch PetscCall(ISSort(valueIS)); 1059566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(valueIS, &Nv)); 1069566063dSJacob Faibussowitsch PetscCall(ISGetIndices(valueIS, &values)); 107adcf2cb4SMatthew G. Knepley /* Construct adaptation label */ 1089566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 1099566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dmCur, 0, &cStart, &cEnd)); 110adcf2cb4SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 111adcf2cb4SMatthew G. Knepley PetscReal volume, centroid[3]; 112adcf2cb4SMatthew G. Knepley PetscInt value, vidx; 113adcf2cb4SMatthew G. Knepley 1149566063dSJacob Faibussowitsch PetscCall(DMPlexComputeCellGeometryFVM(dmCur, c, &volume, centroid, NULL)); 1159566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(label, c, &value)); 116adcf2cb4SMatthew G. Knepley if (value < 0) continue; 1179566063dSJacob Faibussowitsch PetscCall(PetscFindInt(value, Nv, values, &vidx)); 11863a3b9bcSJacob 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); 1199371c9d4SSatish Balay if (volume > volConst[vidx]) { 1209371c9d4SSatish Balay PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE)); 1219371c9d4SSatish Balay ++nAdaptLoc[0]; 1229371c9d4SSatish Balay } 1239371c9d4SSatish Balay if (volume < volConst[vidx] * ratio) { 1249371c9d4SSatish Balay PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_COARSEN)); 1259371c9d4SSatish Balay ++nAdaptLoc[1]; 1269371c9d4SSatish Balay } 127adcf2cb4SMatthew G. Knepley } 1289566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(valueIS, &values)); 1299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&valueIS)); 130462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&nAdaptLoc, &nAdapt, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dmCur))); 131adcf2cb4SMatthew G. Knepley if (nAdapt[0]) { 13263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(dmCur, "Adapted mesh, marking %" PetscInt_FMT " cells for refinement, and %" PetscInt_FMT " cells for coarsening\n", nAdapt[0], nAdapt[1])); 1339566063dSJacob Faibussowitsch PetscCall(DMAdaptLabel(dmCur, adaptLabel, &dmAdapt)); 1349566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmCur)); 1359566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dmAdapt, NULL, "-adapt_dm_view")); 136adcf2cb4SMatthew G. Knepley dmCur = dmAdapt; 137adcf2cb4SMatthew G. Knepley adapt = PETSC_TRUE; 138adcf2cb4SMatthew G. Knepley } 1399566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&adaptLabel)); 140adcf2cb4SMatthew G. Knepley } 1419566063dSJacob Faibussowitsch PetscCall(PetscFree(volConst)); 142adcf2cb4SMatthew G. Knepley *dm = dmCur; 1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 144adcf2cb4SMatthew G. Knepley } 145adcf2cb4SMatthew G. Knepley 146d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 147d71ae5a4SJacob Faibussowitsch { 148adcf2cb4SMatthew G. Knepley PetscInt dim; 149adcf2cb4SMatthew G. Knepley 150adcf2cb4SMatthew G. Knepley PetscFunctionBeginUser; 151adcf2cb4SMatthew G. Knepley /* Create top surface */ 1529566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1539566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1549566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "init_")); 1559566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, NULL)); 157adcf2cb4SMatthew G. Knepley /* Adapt surface */ 1589566063dSJacob Faibussowitsch PetscCall(AdaptMesh(dm, user)); 159adcf2cb4SMatthew G. Knepley /* Extrude surface to get volume mesh */ 1609566063dSJacob Faibussowitsch PetscCall(DMGetDimension(*dm, &dim)); 1619566063dSJacob Faibussowitsch PetscCall(DMLocalizeCoordinates(*dm)); 1629566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, "Mesh")); 1639566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm)); 1649566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 166adcf2cb4SMatthew G. Knepley } 167adcf2cb4SMatthew G. Knepley 168d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 169d71ae5a4SJacob Faibussowitsch { 170adcf2cb4SMatthew G. Knepley DM dm; 171adcf2cb4SMatthew G. Knepley AppCtx user; 172adcf2cb4SMatthew G. Knepley 173327415f7SBarry Smith PetscFunctionBeginUser; 1749566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1759566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 1769566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 1779566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 1789566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 179b122ec5aSJacob Faibussowitsch return 0; 180adcf2cb4SMatthew G. Knepley } 181adcf2cb4SMatthew G. Knepley 182adcf2cb4SMatthew G. Knepley /*TEST 183adcf2cb4SMatthew G. Knepley 184adcf2cb4SMatthew G. Knepley test: 185adcf2cb4SMatthew G. Knepley suffix: 0 186adcf2cb4SMatthew G. Knepley requires: triangle 187d410b0cfSMatthew G. Knepley args: -init_dm_plex_dim 2 -init_dm_plex_box_faces 1,1 -dm_extrude 1 -dm_view 188adcf2cb4SMatthew G. Knepley 189adcf2cb4SMatthew G. Knepley test: # Regularly refine the surface before extrusion 190adcf2cb4SMatthew G. Knepley suffix: 1 191adcf2cb4SMatthew G. Knepley requires: triangle 192d410b0cfSMatthew G. Knepley args: -init_dm_plex_dim 2 -init_dm_refine 2 -dm_extrude 1 -dm_view 193adcf2cb4SMatthew G. Knepley 194adcf2cb4SMatthew G. Knepley test: # Parallel run 195adcf2cb4SMatthew G. Knepley suffix: 2 196adcf2cb4SMatthew G. Knepley requires: triangle 197adcf2cb4SMatthew G. Knepley nsize: 5 198e600fa54SMatthew G. Knepley args: -init_dm_plex_dim 2 -init_dm_refine 3 -petscpartitioner_type simple -dm_extrude 3 -dm_view 199adcf2cb4SMatthew G. Knepley 200adcf2cb4SMatthew G. Knepley test: # adaptively refine the surface before extrusion 201adcf2cb4SMatthew G. Knepley suffix: 3 202adcf2cb4SMatthew G. Knepley requires: triangle 203d410b0cfSMatthew 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 204*e0008caeSPierre Jolivet output_file: output/empty.out 205adcf2cb4SMatthew G. Knepley 206adcf2cb4SMatthew G. Knepley TEST*/ 207