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