1568a717dSMatthew G. Knepley static const char help[] = "Tests for adaptive refinement";
2568a717dSMatthew G. Knepley
3568a717dSMatthew G. Knepley #include <petscdmplex.h>
4*e8b07e64SMatthew G. Knepley #include <petscdmplextransform.h>
5568a717dSMatthew G. Knepley
6568a717dSMatthew G. Knepley typedef struct {
7568a717dSMatthew G. Knepley PetscBool metric; /* Flag to use metric adaptation, instead of tagging */
8568a717dSMatthew G. Knepley PetscInt *refcell; /* A cell to be refined on each process */
9568a717dSMatthew G. Knepley } AppCtx;
10568a717dSMatthew G. Knepley
ProcessOptions(MPI_Comm comm,AppCtx * options)11d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12d71ae5a4SJacob Faibussowitsch {
13568a717dSMatthew G. Knepley PetscMPIInt size;
14568a717dSMatthew G. Knepley PetscInt n;
15568a717dSMatthew G. Knepley
16568a717dSMatthew G. Knepley PetscFunctionBeginUser;
17568a717dSMatthew G. Knepley options->metric = PETSC_FALSE;
189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
199566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(size, &options->refcell));
20568a717dSMatthew G. Knepley n = size;
21568a717dSMatthew G. Knepley
22d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Parallel Mesh Adaptation Options", "DMPLEX");
239566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-metric", "Flag for metric refinement", "ex41.c", options->metric, &options->metric, NULL));
249566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-refcell", "The cell to be refined", "ex41.c", options->refcell, &n, NULL));
2563a3b9bcSJacob Faibussowitsch 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);
26d0609cedSBarry Smith PetscOptionsEnd();
273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
28568a717dSMatthew G. Knepley }
29568a717dSMatthew G. Knepley
CreateMesh(MPI_Comm comm,AppCtx * ctx,DM * dm)30d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
31d71ae5a4SJacob Faibussowitsch {
32568a717dSMatthew G. Knepley PetscFunctionBegin;
339566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
349566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
359566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
369566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
38568a717dSMatthew G. Knepley }
39568a717dSMatthew G. Knepley
CreateAdaptLabel(DM dm,AppCtx * ctx,DMLabel * adaptLabel)40d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateAdaptLabel(DM dm, AppCtx *ctx, DMLabel *adaptLabel)
41d71ae5a4SJacob Faibussowitsch {
42568a717dSMatthew G. Knepley PetscMPIInt rank;
43568a717dSMatthew G. Knepley
44568a717dSMatthew G. Knepley PetscFunctionBegin;
459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
469566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Adaptation Label", adaptLabel));
479566063dSJacob Faibussowitsch if (ctx->refcell[rank] >= 0) PetscCall(DMLabelSetValue(*adaptLabel, ctx->refcell[rank], DM_ADAPT_REFINE));
483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
49568a717dSMatthew G. Knepley }
50568a717dSMatthew G. Knepley
ConstructRefineTree(DM dm)51*e8b07e64SMatthew G. Knepley static PetscErrorCode ConstructRefineTree(DM dm)
52*e8b07e64SMatthew G. Knepley {
53*e8b07e64SMatthew G. Knepley DMPlexTransform tr;
54*e8b07e64SMatthew G. Knepley DM odm;
55*e8b07e64SMatthew G. Knepley PetscInt cStart, cEnd;
56*e8b07e64SMatthew G. Knepley
57*e8b07e64SMatthew G. Knepley PetscFunctionBegin;
58*e8b07e64SMatthew G. Knepley PetscCall(DMPlexGetTransform(dm, &tr));
59*e8b07e64SMatthew G. Knepley if (!tr) PetscFunctionReturn(PETSC_SUCCESS);
60*e8b07e64SMatthew G. Knepley PetscCall(DMPlexTransformGetDM(tr, &odm));
61*e8b07e64SMatthew G. Knepley PetscCall(DMPlexGetHeightStratum(odm, 0, &cStart, &cEnd));
62*e8b07e64SMatthew G. Knepley for (PetscInt c = cStart; c < cEnd; ++c) {
63*e8b07e64SMatthew G. Knepley DMPolytopeType ct;
64*e8b07e64SMatthew G. Knepley DMPolytopeType *rct;
65*e8b07e64SMatthew G. Knepley PetscInt *rsize, *rcone, *rornt;
66*e8b07e64SMatthew G. Knepley PetscInt Nct, dim, pNew = 0;
67*e8b07e64SMatthew G. Knepley
68*e8b07e64SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, "Cell %" PetscInt_FMT " produced new cells", c));
69*e8b07e64SMatthew G. Knepley PetscCall(DMPlexGetCellType(odm, c, &ct));
70*e8b07e64SMatthew G. Knepley dim = DMPolytopeTypeGetDim(ct);
71*e8b07e64SMatthew G. Knepley PetscCall(DMPlexTransformCellTransform(tr, ct, c, NULL, &Nct, &rct, &rsize, &rcone, &rornt));
72*e8b07e64SMatthew G. Knepley for (PetscInt n = 0; n < Nct; ++n) {
73*e8b07e64SMatthew G. Knepley if (DMPolytopeTypeGetDim(rct[n]) != dim) continue;
74*e8b07e64SMatthew G. Knepley for (PetscInt r = 0; r < rsize[n]; ++r) {
75*e8b07e64SMatthew G. Knepley PetscCall(DMPlexTransformGetTargetPoint(tr, ct, rct[n], c, r, &pNew));
76*e8b07e64SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT, pNew));
77*e8b07e64SMatthew G. Knepley }
78*e8b07e64SMatthew G. Knepley }
79*e8b07e64SMatthew G. Knepley PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n"));
80*e8b07e64SMatthew G. Knepley }
81*e8b07e64SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
82*e8b07e64SMatthew G. Knepley }
83*e8b07e64SMatthew G. Knepley
main(int argc,char ** argv)84d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
85d71ae5a4SJacob Faibussowitsch {
86568a717dSMatthew G. Knepley DM dm, dma;
87568a717dSMatthew G. Knepley DMLabel adaptLabel;
88568a717dSMatthew G. Knepley AppCtx ctx;
89568a717dSMatthew G. Knepley
90327415f7SBarry Smith PetscFunctionBeginUser;
919566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
929566063dSJacob Faibussowitsch PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
939566063dSJacob Faibussowitsch PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm));
949566063dSJacob Faibussowitsch PetscCall(CreateAdaptLabel(dm, &ctx, &adaptLabel));
959566063dSJacob Faibussowitsch PetscCall(DMAdaptLabel(dm, adaptLabel, &dma));
969566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dma, "Adapted Mesh"));
979566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&adaptLabel));
989566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
999566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dma, NULL, "-adapt_dm_view"));
100*e8b07e64SMatthew G. Knepley PetscCall(ConstructRefineTree(dma));
1019566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dma));
1029566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx.refcell));
1039566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
104b122ec5aSJacob Faibussowitsch return 0;
105568a717dSMatthew G. Knepley }
106568a717dSMatthew G. Knepley
107568a717dSMatthew G. Knepley /*TEST
108568a717dSMatthew G. Knepley
109012bc364SMatthew G. Knepley testset:
110c0517cd5SMatthew G. Knepley args: -dm_adaptor cellrefiner -dm_plex_transform_type refine_sbr
111012bc364SMatthew G. Knepley
112568a717dSMatthew G. Knepley test:
113568a717dSMatthew G. Knepley suffix: 0
114568a717dSMatthew G. Knepley requires: triangle
115012bc364SMatthew G. Knepley args: -dm_view -adapt_dm_view
116568a717dSMatthew G. Knepley
117568a717dSMatthew G. Knepley test:
118568a717dSMatthew G. Knepley suffix: 1
119568a717dSMatthew G. Knepley requires: triangle
120012bc364SMatthew G. Knepley args: -dm_coord_space 0 -refcell 2 -dm_view ::ascii_info_detail -adapt_dm_view ::ascii_info_detail
121568a717dSMatthew G. Knepley
122568a717dSMatthew G. Knepley test:
123*e8b07e64SMatthew G. Knepley suffix: 1_save
124*e8b07e64SMatthew G. Knepley requires: triangle
125*e8b07e64SMatthew G. Knepley args: -refcell 2 -dm_plex_save_transform -dm_view -adapt_dm_view
126*e8b07e64SMatthew G. Knepley
127*e8b07e64SMatthew G. Knepley test:
128568a717dSMatthew G. Knepley suffix: 2
129568a717dSMatthew G. Knepley requires: triangle
130568a717dSMatthew G. Knepley nsize: 2
131e600fa54SMatthew G. Knepley args: -refcell 2,-1 -petscpartitioner_type simple -dm_view -adapt_dm_view
132568a717dSMatthew G. Knepley
133568a717dSMatthew G. Knepley TEST*/
134