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 PetscErrorCode ierr; 15 16 PetscFunctionBeginUser; 17 options->metric = PETSC_FALSE; 18 ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 19 ierr = PetscCalloc1(size, &options->refcell);CHKERRQ(ierr); 20 n = size; 21 22 ierr = PetscOptionsBegin(comm, "", "Parallel Mesh Adaptation Options", "DMPLEX");CHKERRQ(ierr); 23 ierr = PetscOptionsBool("-metric", "Flag for metric refinement", "ex41.c", options->metric, &options->metric, NULL);CHKERRQ(ierr); 24 ierr = PetscOptionsIntArray("-refcell", "The cell to be refined", "ex41.c", options->refcell, &n, NULL);CHKERRQ(ierr); 25 if (n && n != size) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Only gave %D cells to refine, must give one for all %D processes", n, size); 26 ierr = PetscOptionsEnd();CHKERRQ(ierr); 27 PetscFunctionReturn(0); 28 } 29 30 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 31 { 32 PetscErrorCode ierr; 33 34 PetscFunctionBegin; 35 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 36 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 37 ierr = DMSetFromOptions(*dm);CHKERRQ(ierr); 38 ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr); 39 PetscFunctionReturn(0); 40 } 41 42 static PetscErrorCode CreateAdaptLabel(DM dm, AppCtx *ctx, DMLabel *adaptLabel) 43 { 44 PetscMPIInt rank; 45 PetscErrorCode ierr; 46 47 PetscFunctionBegin; 48 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRMPI(ierr); 49 ierr = DMLabelCreate(PETSC_COMM_SELF, "Adaptation Label", adaptLabel);CHKERRQ(ierr); 50 if (ctx->refcell[rank] >= 0) {ierr = DMLabelSetValue(*adaptLabel, ctx->refcell[rank], DM_ADAPT_REFINE);CHKERRQ(ierr);} 51 PetscFunctionReturn(0); 52 } 53 54 int main(int argc, char **argv) 55 { 56 DM dm, dma; 57 DMLabel adaptLabel; 58 AppCtx ctx; 59 PetscErrorCode ierr; 60 61 ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr; 62 ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr); 63 ierr = CreateMesh(PETSC_COMM_WORLD, &ctx, &dm);CHKERRQ(ierr); 64 ierr = CreateAdaptLabel(dm, &ctx, &adaptLabel);CHKERRQ(ierr); 65 ierr = DMAdaptLabel(dm, adaptLabel, &dma);CHKERRQ(ierr); 66 ierr = PetscObjectSetName((PetscObject) dma, "Adapted Mesh");CHKERRQ(ierr); 67 ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr); 68 ierr = DMDestroy(&dm);CHKERRQ(ierr); 69 ierr = DMViewFromOptions(dma, NULL, "-adapt_dm_view");CHKERRQ(ierr); 70 ierr = DMDestroy(&dma);CHKERRQ(ierr); 71 ierr = PetscFree(ctx.refcell);CHKERRQ(ierr); 72 ierr = PetscFinalize(); 73 return ierr; 74 } 75 76 /*TEST 77 78 testset: 79 args: -dm_plex_adaptor cellrefiner -dm_plex_transform_type refine_sbr 80 81 test: 82 suffix: 0 83 requires: triangle 84 args: -dm_view -adapt_dm_view 85 86 test: 87 suffix: 1 88 requires: triangle 89 args: -dm_coord_space 0 -refcell 2 -dm_view ::ascii_info_detail -adapt_dm_view ::ascii_info_detail 90 91 test: 92 suffix: 2 93 requires: triangle 94 nsize: 2 95 args: -refcell 2,-1 -dm_distribute -petscpartitioner_type simple -dm_view -adapt_dm_view 96 97 TEST*/ 98