1503c0ea9SMatthew G. Knepley static const char help[] = "Tests for mesh extrusion"; 2503c0ea9SMatthew G. Knepley 3503c0ea9SMatthew G. Knepley #include <petscdmplex.h> 4503c0ea9SMatthew G. Knepley 5503c0ea9SMatthew G. Knepley typedef struct { 6503c0ea9SMatthew G. Knepley char bdLabel[PETSC_MAX_PATH_LEN]; /* The boundary label name */ 7503c0ea9SMatthew G. Knepley PetscInt Nbd; /* The number of boundary markers to extrude, 0 for all */ 8503c0ea9SMatthew G. Knepley PetscInt bd[64]; /* The boundary markers to be extruded */ 9503c0ea9SMatthew G. Knepley } AppCtx; 10503c0ea9SMatthew G. Knepley 11503c0ea9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode pyramidNormal(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 12503c0ea9SMatthew G. Knepley 13503c0ea9SMatthew G. Knepley /* The pyramid apex is at (0.5, 0.5, -1) */ 14503c0ea9SMatthew G. Knepley PetscErrorCode pyramidNormal(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt r, PetscScalar u[], void *ctx) 15503c0ea9SMatthew G. Knepley { 16503c0ea9SMatthew G. Knepley PetscReal apex[3] = {0.5, 0.5, -1.0}; 17503c0ea9SMatthew G. Knepley PetscInt d; 18503c0ea9SMatthew G. Knepley 19503c0ea9SMatthew G. Knepley for (d = 0; d < dim; ++d) u[d] = x[d] - apex[d]; 20503c0ea9SMatthew G. Knepley for (d = dim; d < 3; ++d) u[d] = 0.0 - apex[d]; 21503c0ea9SMatthew G. Knepley return 0; 22503c0ea9SMatthew G. Knepley } 23503c0ea9SMatthew G. Knepley 24503c0ea9SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 25503c0ea9SMatthew G. Knepley { 26503c0ea9SMatthew G. Knepley PetscInt n = 64; 27503c0ea9SMatthew G. Knepley PetscBool flg; 28503c0ea9SMatthew G. Knepley 29503c0ea9SMatthew G. Knepley PetscFunctionBeginUser; 30503c0ea9SMatthew G. Knepley PetscCall(PetscStrcpy(options->bdLabel, "marker")); 31*d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Parallel Mesh Adaptation Options", "DMPLEX"); 32503c0ea9SMatthew G. Knepley PetscCall(PetscOptionsString("-label", "The boundary label name", "ex44.c", options->bdLabel, options->bdLabel, sizeof(options->bdLabel), NULL)); 33503c0ea9SMatthew G. Knepley PetscCall(PetscOptionsIntArray("-bd", "The boundaries to be extruded", "ex44.c", options->bd, &n, &flg)); 34503c0ea9SMatthew G. Knepley options->Nbd = flg ? n : 0; 35*d0609cedSBarry Smith PetscOptionsEnd(); 36503c0ea9SMatthew G. Knepley PetscFunctionReturn(0); 37503c0ea9SMatthew G. Knepley } 38503c0ea9SMatthew G. Knepley 39503c0ea9SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 40503c0ea9SMatthew G. Knepley { 41503c0ea9SMatthew G. Knepley PetscFunctionBegin; 42503c0ea9SMatthew G. Knepley PetscCall(DMCreate(comm, dm)); 43503c0ea9SMatthew G. Knepley PetscCall(DMSetType(*dm, DMPLEX)); 44503c0ea9SMatthew G. Knepley PetscCall(DMSetFromOptions(*dm)); 45503c0ea9SMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 46503c0ea9SMatthew G. Knepley PetscFunctionReturn(0); 47503c0ea9SMatthew G. Knepley } 48503c0ea9SMatthew G. Knepley 49503c0ea9SMatthew G. Knepley static PetscErrorCode CreateAdaptLabel(DM dm, AppCtx *ctx, DMLabel *adaptLabel) 50503c0ea9SMatthew G. Knepley { 51503c0ea9SMatthew G. Knepley DMLabel label; 52503c0ea9SMatthew G. Knepley PetscInt b; 53503c0ea9SMatthew G. Knepley 54503c0ea9SMatthew G. Knepley PetscFunctionBegin; 55503c0ea9SMatthew G. Knepley if (!ctx->Nbd) {*adaptLabel = NULL; PetscFunctionReturn(0);} 56503c0ea9SMatthew G. Knepley PetscCall(DMGetLabel(dm, ctx->bdLabel, &label)); 57503c0ea9SMatthew G. Knepley PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Adaptation Label", adaptLabel)); 58503c0ea9SMatthew G. Knepley for (b = 0; b < ctx->Nbd; ++b) { 59503c0ea9SMatthew G. Knepley IS bdIS; 60503c0ea9SMatthew G. Knepley const PetscInt *points; 61503c0ea9SMatthew G. Knepley PetscInt n, i; 62503c0ea9SMatthew G. Knepley 63503c0ea9SMatthew G. Knepley PetscCall(DMLabelGetStratumIS(label, ctx->bd[b], &bdIS)); 64503c0ea9SMatthew G. Knepley if (!bdIS) continue; 65503c0ea9SMatthew G. Knepley PetscCall(ISGetLocalSize(bdIS, &n)); 66503c0ea9SMatthew G. Knepley PetscCall(ISGetIndices(bdIS, &points)); 67503c0ea9SMatthew G. Knepley for (i = 0; i < n; ++i) {PetscCall(DMLabelSetValue(*adaptLabel, points[i], DM_ADAPT_REFINE));} 68503c0ea9SMatthew G. Knepley PetscCall(ISRestoreIndices(bdIS, &points)); 69503c0ea9SMatthew G. Knepley PetscCall(ISDestroy(&bdIS)); 70503c0ea9SMatthew G. Knepley } 71503c0ea9SMatthew G. Knepley PetscFunctionReturn(0); 72503c0ea9SMatthew G. Knepley } 73503c0ea9SMatthew G. Knepley 74503c0ea9SMatthew G. Knepley int main(int argc, char **argv) 75503c0ea9SMatthew G. Knepley { 76503c0ea9SMatthew G. Knepley DM dm, dma; 77503c0ea9SMatthew G. Knepley DMLabel adaptLabel; 78503c0ea9SMatthew G. Knepley AppCtx ctx; 79503c0ea9SMatthew G. Knepley 80503c0ea9SMatthew G. Knepley PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 81503c0ea9SMatthew G. Knepley PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 82503c0ea9SMatthew G. Knepley PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm)); 83503c0ea9SMatthew G. Knepley PetscCall(CreateAdaptLabel(dm, &ctx, &adaptLabel)); 84503c0ea9SMatthew G. Knepley if (adaptLabel) {PetscCall(DMAdaptLabel(dm, adaptLabel, &dma));} 85503c0ea9SMatthew G. Knepley else {PetscCall(DMExtrude(dm, 3, &dma));} 86503c0ea9SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject) dma, "Adapted Mesh")); 87503c0ea9SMatthew G. Knepley PetscCall(DMLabelDestroy(&adaptLabel)); 88503c0ea9SMatthew G. Knepley PetscCall(DMDestroy(&dm)); 89503c0ea9SMatthew G. Knepley PetscCall(DMViewFromOptions(dma, NULL, "-adapt_dm_view")); 90503c0ea9SMatthew G. Knepley PetscCall(DMDestroy(&dma)); 91503c0ea9SMatthew G. Knepley PetscCall(PetscFinalize()); 92503c0ea9SMatthew G. Knepley return 0; 93503c0ea9SMatthew G. Knepley } 94503c0ea9SMatthew G. Knepley 95503c0ea9SMatthew G. Knepley /*TEST 96503c0ea9SMatthew G. Knepley 97503c0ea9SMatthew G. Knepley test: 98503c0ea9SMatthew G. Knepley suffix: tri_tensor_0 99503c0ea9SMatthew G. Knepley requires: triangle 100503c0ea9SMatthew G. Knepley args: -dm_plex_transform_extrude_use_tensor {{0 1}separate output} \ 101503c0ea9SMatthew G. Knepley -dm_view -adapt_dm_view -dm_plex_check_all 102503c0ea9SMatthew G. Knepley 103503c0ea9SMatthew G. Knepley test: 104503c0ea9SMatthew G. Knepley suffix: quad_tensor_0 105503c0ea9SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_transform_extrude_use_tensor {{0 1}separate output} \ 106503c0ea9SMatthew G. Knepley -dm_view -adapt_dm_view -dm_plex_check_all 107503c0ea9SMatthew G. Knepley 108503c0ea9SMatthew G. Knepley test: 109503c0ea9SMatthew G. Knepley suffix: quad_normal_0 110503c0ea9SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_transform_extrude_normal 0,1,1 \ 111503c0ea9SMatthew G. Knepley -dm_view -adapt_dm_view -dm_plex_check_all 112503c0ea9SMatthew G. Knepley 113503c0ea9SMatthew G. Knepley test: 114503c0ea9SMatthew G. Knepley suffix: quad_normal_1 115503c0ea9SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_transform_extrude_normal_function pyramidNormal \ 116503c0ea9SMatthew G. Knepley -dm_view -adapt_dm_view -dm_plex_check_all 117503c0ea9SMatthew G. Knepley 118503c0ea9SMatthew G. Knepley test: 119503c0ea9SMatthew G. Knepley suffix: quad_symmetric_0 120503c0ea9SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_transform_extrude_symmetric \ 121503c0ea9SMatthew G. Knepley -dm_view -adapt_dm_view -dm_plex_check_all 122503c0ea9SMatthew G. Knepley 123503c0ea9SMatthew G. Knepley testset: 124503c0ea9SMatthew G. Knepley args: -dm_adaptor cellrefiner -dm_plex_transform_type extrude \ 125503c0ea9SMatthew G. Knepley -dm_view -adapt_dm_view 126503c0ea9SMatthew G. Knepley 127503c0ea9SMatthew G. Knepley test: 128503c0ea9SMatthew G. Knepley suffix: quad_adapt_0 129503c0ea9SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -dm_plex_separate_marker -bd 1,3 130503c0ea9SMatthew G. Knepley 131503c0ea9SMatthew G. Knepley TEST*/ 132