1 static const char help[] = "Tests for adaptive refinement"; 2 3 #include <petscdmplex.h> 4 5 typedef struct { 6 PetscBool metric; /* Flag to use metric adaptation, instead of tagging */ 7 PetscInt *refcell; /* A cell to be refined on each process */ 8 } AppCtx; 9 10 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 11 { 12 PetscMPIInt size; 13 PetscInt n; 14 15 PetscFunctionBeginUser; 16 options->metric = PETSC_FALSE; 17 PetscCallMPI(MPI_Comm_size(comm, &size)); 18 PetscCall(PetscCalloc1(size, &options->refcell)); 19 n = size; 20 21 PetscOptionsBegin(comm, "", "Parallel Mesh Adaptation Options", "DMPLEX"); 22 PetscCall(PetscOptionsBool("-metric", "Flag for metric refinement", "ex41.c", options->metric, &options->metric, NULL)); 23 PetscCall(PetscOptionsIntArray("-refcell", "The cell to be refined", "ex41.c", options->refcell, &n, NULL)); 24 if (n) PetscCheck(n == size, comm, PETSC_ERR_ARG_SIZ, "Only gave %" PetscInt_FMT " cells to refine, must give one for all %d processes", n, size); 25 PetscOptionsEnd(); 26 PetscFunctionReturn(0); 27 } 28 29 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 30 { 31 PetscFunctionBegin; 32 PetscCall(DMCreate(comm, dm)); 33 PetscCall(DMSetType(*dm, DMPLEX)); 34 PetscCall(DMSetFromOptions(*dm)); 35 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 36 PetscFunctionReturn(0); 37 } 38 39 static PetscErrorCode CreateAdaptLabel(DM dm, AppCtx *ctx, DMLabel *adaptLabel) 40 { 41 PetscMPIInt rank; 42 43 PetscFunctionBegin; 44 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 45 PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Adaptation Label", adaptLabel)); 46 if (ctx->refcell[rank] >= 0) PetscCall(DMLabelSetValue(*adaptLabel, ctx->refcell[rank], DM_ADAPT_REFINE)); 47 PetscFunctionReturn(0); 48 } 49 50 int main(int argc, char **argv) 51 { 52 DM dm, dma; 53 DMLabel adaptLabel; 54 AppCtx ctx; 55 56 PetscFunctionBeginUser; 57 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 58 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 59 PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm)); 60 PetscCall(CreateAdaptLabel(dm, &ctx, &adaptLabel)); 61 PetscCall(DMAdaptLabel(dm, adaptLabel, &dma)); 62 PetscCall(PetscObjectSetName((PetscObject)dma, "Adapted Mesh")); 63 PetscCall(DMLabelDestroy(&adaptLabel)); 64 PetscCall(DMDestroy(&dm)); 65 PetscCall(DMViewFromOptions(dma, NULL, "-adapt_dm_view")); 66 PetscCall(DMDestroy(&dma)); 67 PetscCall(PetscFree(ctx.refcell)); 68 PetscCall(PetscFinalize()); 69 return 0; 70 } 71 72 /*TEST 73 74 testset: 75 args: -dm_adaptor cellrefiner -dm_plex_transform_type refine_sbr 76 77 test: 78 suffix: 0 79 requires: triangle 80 args: -dm_view -adapt_dm_view 81 82 test: 83 suffix: 1 84 requires: triangle 85 args: -dm_coord_space 0 -refcell 2 -dm_view ::ascii_info_detail -adapt_dm_view ::ascii_info_detail 86 87 test: 88 suffix: 2 89 requires: triangle 90 nsize: 2 91 args: -refcell 2,-1 -petscpartitioner_type simple -dm_view -adapt_dm_view 92 93 TEST*/ 94