1568a717dSMatthew G. Knepley static const char help[] = "Tests for adaptive refinement"; 2568a717dSMatthew G. Knepley 3568a717dSMatthew G. Knepley #include <petscdmplex.h> 4568a717dSMatthew G. Knepley 5568a717dSMatthew G. Knepley typedef struct { 6568a717dSMatthew G. Knepley PetscBool metric; /* Flag to use metric adaptation, instead of tagging */ 7568a717dSMatthew G. Knepley PetscInt *refcell; /* A cell to be refined on each process */ 8568a717dSMatthew G. Knepley } AppCtx; 9568a717dSMatthew G. Knepley 10568a717dSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 11568a717dSMatthew G. Knepley { 12568a717dSMatthew G. Knepley PetscMPIInt size; 13568a717dSMatthew G. Knepley PetscInt n; 14568a717dSMatthew G. Knepley PetscErrorCode ierr; 15568a717dSMatthew G. Knepley 16568a717dSMatthew G. Knepley PetscFunctionBeginUser; 17568a717dSMatthew G. Knepley options->metric = PETSC_FALSE; 18*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm, &size)); 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(size, &options->refcell)); 20568a717dSMatthew G. Knepley n = size; 21568a717dSMatthew G. Knepley 22568a717dSMatthew G. Knepley ierr = PetscOptionsBegin(comm, "", "Parallel Mesh Adaptation Options", "DMPLEX");CHKERRQ(ierr); 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-metric", "Flag for metric refinement", "ex41.c", options->metric, &options->metric, NULL)); 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsIntArray("-refcell", "The cell to be refined", "ex41.c", options->refcell, &n, NULL)); 252c71b3e2SJacob Faibussowitsch PetscCheckFalse(n && n != size,comm, PETSC_ERR_ARG_SIZ, "Only gave %D cells to refine, must give one for all %D processes", n, size); 261e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 27568a717dSMatthew G. Knepley PetscFunctionReturn(0); 28568a717dSMatthew G. Knepley } 29568a717dSMatthew G. Knepley 30568a717dSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 31568a717dSMatthew G. Knepley { 32568a717dSMatthew G. Knepley PetscFunctionBegin; 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(*dm)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 37568a717dSMatthew G. Knepley PetscFunctionReturn(0); 38568a717dSMatthew G. Knepley } 39568a717dSMatthew G. Knepley 40568a717dSMatthew G. Knepley static PetscErrorCode CreateAdaptLabel(DM dm, AppCtx *ctx, DMLabel *adaptLabel) 41568a717dSMatthew G. Knepley { 42568a717dSMatthew G. Knepley PetscMPIInt rank; 43568a717dSMatthew G. Knepley 44568a717dSMatthew G. Knepley PetscFunctionBegin; 45*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "Adaptation Label", adaptLabel)); 47*5f80ce2aSJacob Faibussowitsch if (ctx->refcell[rank] >= 0) CHKERRQ(DMLabelSetValue(*adaptLabel, ctx->refcell[rank], DM_ADAPT_REFINE)); 48568a717dSMatthew G. Knepley PetscFunctionReturn(0); 49568a717dSMatthew G. Knepley } 50568a717dSMatthew G. Knepley 51568a717dSMatthew G. Knepley int main(int argc, char **argv) 52568a717dSMatthew G. Knepley { 53568a717dSMatthew G. Knepley DM dm, dma; 54568a717dSMatthew G. Knepley DMLabel adaptLabel; 55568a717dSMatthew G. Knepley AppCtx ctx; 56568a717dSMatthew G. Knepley PetscErrorCode ierr; 57568a717dSMatthew G. Knepley 58568a717dSMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr; 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(CreateAdaptLabel(dm, &ctx, &adaptLabel)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMAdaptLabel(dm, adaptLabel, &dma)); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) dma, "Adapted Mesh")); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelDestroy(&adaptLabel)); 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dm)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMViewFromOptions(dma, NULL, "-adapt_dm_view")); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&dma)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ctx.refcell)); 69568a717dSMatthew G. Knepley ierr = PetscFinalize(); 70568a717dSMatthew G. Knepley return ierr; 71568a717dSMatthew G. Knepley } 72568a717dSMatthew G. Knepley 73568a717dSMatthew G. Knepley /*TEST 74568a717dSMatthew G. Knepley 75012bc364SMatthew G. Knepley testset: 76c0517cd5SMatthew G. Knepley args: -dm_adaptor cellrefiner -dm_plex_transform_type refine_sbr 77012bc364SMatthew G. Knepley 78568a717dSMatthew G. Knepley test: 79568a717dSMatthew G. Knepley suffix: 0 80568a717dSMatthew G. Knepley requires: triangle 81012bc364SMatthew G. Knepley args: -dm_view -adapt_dm_view 82568a717dSMatthew G. Knepley 83568a717dSMatthew G. Knepley test: 84568a717dSMatthew G. Knepley suffix: 1 85568a717dSMatthew G. Knepley requires: triangle 86012bc364SMatthew G. Knepley args: -dm_coord_space 0 -refcell 2 -dm_view ::ascii_info_detail -adapt_dm_view ::ascii_info_detail 87568a717dSMatthew G. Knepley 88568a717dSMatthew G. Knepley test: 89568a717dSMatthew G. Knepley suffix: 2 90568a717dSMatthew G. Knepley requires: triangle 91568a717dSMatthew G. Knepley nsize: 2 92e600fa54SMatthew G. Knepley args: -refcell 2,-1 -petscpartitioner_type simple -dm_view -adapt_dm_view 93568a717dSMatthew G. Knepley 94568a717dSMatthew G. Knepley TEST*/ 95