xref: /petsc/src/dm/impls/plex/tests/ex44.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1 static const char help[] = "Tests for mesh extrusion";
2 
3 #include <petscdmplex.h>
4 
5 typedef struct {
6   char     bdLabel[PETSC_MAX_PATH_LEN]; /* The boundary label name */
7   PetscInt Nbd;                         /* The number of boundary markers to extrude, 0 for all */
8   PetscInt bd[64];                      /* The boundary markers to be extruded */
9 } AppCtx;
10 
11 PETSC_EXTERN PetscErrorCode pyramidNormal(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *);
12 
13 /* The pyramid apex is at (0.5, 0.5, -1) */
14 PetscErrorCode pyramidNormal(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt r, PetscScalar u[], void *ctx) {
15   PetscReal apex[3] = {0.5, 0.5, -1.0};
16   PetscInt  d;
17 
18   for (d = 0; d < dim; ++d) u[d] = x[d] - apex[d];
19   for (d = dim; d < 3; ++d) u[d] = 0.0 - apex[d];
20   return 0;
21 }
22 
23 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
24   PetscInt  n = 64;
25   PetscBool flg;
26 
27   PetscFunctionBeginUser;
28   PetscCall(PetscStrcpy(options->bdLabel, "marker"));
29   PetscOptionsBegin(comm, "", "Parallel Mesh Adaptation Options", "DMPLEX");
30   PetscCall(PetscOptionsString("-label", "The boundary label name", "ex44.c", options->bdLabel, options->bdLabel, sizeof(options->bdLabel), NULL));
31   PetscCall(PetscOptionsIntArray("-bd", "The boundaries to be extruded", "ex44.c", options->bd, &n, &flg));
32   options->Nbd = flg ? n : 0;
33   PetscOptionsEnd();
34   PetscFunctionReturn(0);
35 }
36 
37 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) {
38   PetscFunctionBegin;
39   PetscCall(DMCreate(comm, dm));
40   PetscCall(DMSetType(*dm, DMPLEX));
41   PetscCall(DMSetFromOptions(*dm));
42   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
43   PetscFunctionReturn(0);
44 }
45 
46 static PetscErrorCode CreateAdaptLabel(DM dm, AppCtx *ctx, DMLabel *adaptLabel) {
47   DMLabel  label;
48   PetscInt b;
49 
50   PetscFunctionBegin;
51   if (!ctx->Nbd) {
52     *adaptLabel = NULL;
53     PetscFunctionReturn(0);
54   }
55   PetscCall(DMGetLabel(dm, ctx->bdLabel, &label));
56   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Adaptation Label", adaptLabel));
57   for (b = 0; b < ctx->Nbd; ++b) {
58     IS              bdIS;
59     const PetscInt *points;
60     PetscInt        n, i;
61 
62     PetscCall(DMLabelGetStratumIS(label, ctx->bd[b], &bdIS));
63     if (!bdIS) continue;
64     PetscCall(ISGetLocalSize(bdIS, &n));
65     PetscCall(ISGetIndices(bdIS, &points));
66     for (i = 0; i < n; ++i) PetscCall(DMLabelSetValue(*adaptLabel, points[i], DM_ADAPT_REFINE));
67     PetscCall(ISRestoreIndices(bdIS, &points));
68     PetscCall(ISDestroy(&bdIS));
69   }
70   PetscFunctionReturn(0);
71 }
72 
73 int main(int argc, char **argv) {
74   DM      dm, dma;
75   DMLabel adaptLabel;
76   AppCtx  ctx;
77 
78   PetscFunctionBeginUser;
79   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
80   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
81   PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm));
82   PetscCall(CreateAdaptLabel(dm, &ctx, &adaptLabel));
83   if (adaptLabel) {
84     PetscCall(DMAdaptLabel(dm, adaptLabel, &dma));
85   } else {
86     PetscCall(DMExtrude(dm, 3, &dma));
87   }
88   PetscCall(PetscObjectSetName((PetscObject)dma, "Adapted Mesh"));
89   PetscCall(DMLabelDestroy(&adaptLabel));
90   PetscCall(DMDestroy(&dm));
91   PetscCall(DMViewFromOptions(dma, NULL, "-adapt_dm_view"));
92   PetscCall(DMDestroy(&dma));
93   PetscCall(PetscFinalize());
94   return 0;
95 }
96 
97 /*TEST
98 
99   test:
100     suffix: tri_tensor_0
101     requires: triangle
102     args: -dm_plex_transform_extrude_use_tensor {{0 1}separate output} \
103           -dm_view -adapt_dm_view -dm_plex_check_all
104 
105   test:
106     suffix: quad_tensor_0
107     args: -dm_plex_simplex 0 -dm_plex_transform_extrude_use_tensor {{0 1}separate output} \
108           -dm_view -adapt_dm_view -dm_plex_check_all
109 
110   test:
111     suffix: quad_normal_0
112     args: -dm_plex_simplex 0 -dm_plex_transform_extrude_normal 0,1,1 \
113           -dm_view -adapt_dm_view -dm_plex_check_all
114 
115   test:
116     suffix: quad_normal_1
117     args: -dm_plex_simplex 0 -dm_plex_transform_extrude_normal_function pyramidNormal \
118           -dm_view -adapt_dm_view -dm_plex_check_all
119 
120   test:
121     suffix: quad_symmetric_0
122     args: -dm_plex_simplex 0 -dm_plex_transform_extrude_symmetric \
123           -dm_view -adapt_dm_view -dm_plex_check_all
124 
125   testset:
126     args: -dm_adaptor cellrefiner -dm_plex_transform_type extrude \
127           -dm_view -adapt_dm_view
128 
129     test:
130       suffix: quad_adapt_0
131       args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -dm_plex_separate_marker -bd 1,3 \
132             -dm_plex_transform_extrude_thickness 0.5
133 
134     test:
135       suffix: tet_adapt_0
136       requires: ctetgen
137       args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_separate_marker -bd 1,3 \
138             -dm_plex_transform_extrude_thickness 0.5
139 
140     test:
141       suffix: hex_adapt_0
142       args: -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_separate_marker -bd 1,3 \
143             -dm_plex_transform_extrude_thickness 0.5
144 
145 TEST*/
146