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