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