xref: /petsc/src/dm/impls/plex/tests/ex41.c (revision ccfb0f9f40a0131988d7995ed9679700dae2a75a)
1 static const char help[] = "Tests for adaptive refinement";
2 
3 #include <petscdmplex.h>
4 #include <petscdmplextransform.h>
5 
6 typedef struct {
7   PetscBool metric;  /* Flag to use metric adaptation, instead of tagging */
8   PetscInt *refcell; /* A cell to be refined on each process */
9 } AppCtx;
10 
11 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12 {
13   PetscMPIInt size;
14   PetscInt    n;
15 
16   PetscFunctionBeginUser;
17   options->metric = PETSC_FALSE;
18   PetscCallMPI(MPI_Comm_size(comm, &size));
19   PetscCall(PetscCalloc1(size, &options->refcell));
20   n = size;
21 
22   PetscOptionsBegin(comm, "", "Parallel Mesh Adaptation Options", "DMPLEX");
23   PetscCall(PetscOptionsBool("-metric", "Flag for metric refinement", "ex41.c", options->metric, &options->metric, NULL));
24   PetscCall(PetscOptionsIntArray("-refcell", "The cell to be refined", "ex41.c", options->refcell, &n, NULL));
25   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);
26   PetscOptionsEnd();
27   PetscFunctionReturn(PETSC_SUCCESS);
28 }
29 
30 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
31 {
32   PetscFunctionBegin;
33   PetscCall(DMCreate(comm, dm));
34   PetscCall(DMSetType(*dm, DMPLEX));
35   PetscCall(DMSetFromOptions(*dm));
36   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
37   PetscFunctionReturn(PETSC_SUCCESS);
38 }
39 
40 static PetscErrorCode CreateAdaptLabel(DM dm, AppCtx *ctx, DMLabel *adaptLabel)
41 {
42   PetscMPIInt rank;
43 
44   PetscFunctionBegin;
45   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
46   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Adaptation Label", adaptLabel));
47   if (ctx->refcell[rank] >= 0) PetscCall(DMLabelSetValue(*adaptLabel, ctx->refcell[rank], DM_ADAPT_REFINE));
48   PetscFunctionReturn(PETSC_SUCCESS);
49 }
50 
51 static PetscErrorCode ConstructRefineTree(DM dm)
52 {
53   DMPlexTransform tr;
54   DM              odm;
55   PetscInt        cStart, cEnd;
56 
57   PetscFunctionBegin;
58   PetscCall(DMPlexGetTransform(dm, &tr));
59   if (!tr) PetscFunctionReturn(PETSC_SUCCESS);
60   PetscCall(DMPlexTransformGetDM(tr, &odm));
61   PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
62   for (PetscInt c = cStart; c < cEnd; ++c) {
63     DMPolytopeType  ct;
64     DMPolytopeType *rct;
65     PetscInt       *rsize, *rcone, *rornt;
66     PetscInt        Nct, dim, pNew = 0;
67 
68     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " produced new cells", c));
69     PetscCall(DMPlexGetCellType(odm, c, &ct));
70     dim = DMPolytopeTypeGetDim(ct);
71     PetscCall(DMPlexTransformCellTransform(tr, ct, c, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
72     for (PetscInt n = 0; n < Nct; ++n) {
73       if (DMPolytopeTypeGetDim(rct[n]) != dim) continue;
74       for (PetscInt r = 0; r < rsize[n]; ++r) {
75         PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[n], c, r, &pNew));
76         PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT, pNew));
77       }
78     }
79     PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
80   }
81   PetscFunctionReturn(PETSC_SUCCESS);
82 }
83 
84 int main(int argc, char **argv)
85 {
86   DM      dm, dma;
87   DMLabel adaptLabel;
88   AppCtx  ctx;
89 
90   PetscFunctionBeginUser;
91   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
92   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
93   PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm));
94   PetscCall(CreateAdaptLabel(dm, &ctx, &adaptLabel));
95   PetscCall(DMAdaptLabel(dm, adaptLabel, &dma));
96   PetscCall(PetscObjectSetName((PetscObject)dma, "Adapted Mesh"));
97   PetscCall(DMLabelDestroy(&adaptLabel));
98   PetscCall(DMDestroy(&dm));
99   PetscCall(DMViewFromOptions(dma, NULL, "-adapt_dm_view"));
100   PetscCall(ConstructRefineTree(dma));
101   PetscCall(DMDestroy(&dma));
102   PetscCall(PetscFree(ctx.refcell));
103   PetscCall(PetscFinalize());
104   return 0;
105 }
106 
107 /*TEST
108 
109   testset:
110     args: -dm_adaptor cellrefiner -dm_plex_transform_type refine_sbr
111 
112     test:
113       suffix: 0
114       requires: triangle
115       args: -dm_view -adapt_dm_view
116 
117     test:
118       suffix: 1
119       requires: triangle
120       args: -dm_coord_space 0 -refcell 2 -dm_view ::ascii_info_detail -adapt_dm_view ::ascii_info_detail
121 
122     test:
123       suffix: 1_save
124       requires: triangle
125       args: -refcell 2 -dm_plex_save_transform -dm_view -adapt_dm_view
126 
127     test:
128       suffix: 2
129       requires: triangle
130       nsize: 2
131       args: -refcell 2,-1 -petscpartitioner_type simple -dm_view -adapt_dm_view
132 
133 TEST*/
134