xref: /petsc/src/dm/impls/plex/tests/ex41.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1568a717dSMatthew G. Knepley static const char help[] = "Tests for adaptive refinement";
2568a717dSMatthew G. Knepley 
3568a717dSMatthew G. Knepley #include <petscdmplex.h>
4568a717dSMatthew G. Knepley 
5568a717dSMatthew G. Knepley typedef struct {
6568a717dSMatthew G. Knepley   PetscBool metric;  /* Flag to use metric adaptation, instead of tagging */
7568a717dSMatthew G. Knepley   PetscInt *refcell; /* A cell to be refined on each process */
8568a717dSMatthew G. Knepley } AppCtx;
9568a717dSMatthew G. Knepley 
10568a717dSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
11568a717dSMatthew G. Knepley {
12568a717dSMatthew G. Knepley   PetscMPIInt size;
13568a717dSMatthew G. Knepley   PetscInt    n;
14568a717dSMatthew G. Knepley 
15568a717dSMatthew G. Knepley   PetscFunctionBeginUser;
16568a717dSMatthew G. Knepley   options->metric  = PETSC_FALSE;
179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
189566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(size, &options->refcell));
19568a717dSMatthew G. Knepley   n = size;
20568a717dSMatthew G. Knepley 
21d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Parallel Mesh Adaptation Options", "DMPLEX");
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-metric", "Flag for metric refinement", "ex41.c", options->metric, &options->metric, NULL));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-refcell", "The cell to be refined", "ex41.c", options->refcell, &n, NULL));
2463a3b9bcSJacob 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);
25d0609cedSBarry Smith   PetscOptionsEnd();
26568a717dSMatthew G. Knepley   PetscFunctionReturn(0);
27568a717dSMatthew G. Knepley }
28568a717dSMatthew G. Knepley 
29568a717dSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
30568a717dSMatthew G. Knepley {
31568a717dSMatthew G. Knepley   PetscFunctionBegin;
329566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
339566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
349566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
359566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
36568a717dSMatthew G. Knepley   PetscFunctionReturn(0);
37568a717dSMatthew G. Knepley }
38568a717dSMatthew G. Knepley 
39568a717dSMatthew G. Knepley static PetscErrorCode CreateAdaptLabel(DM dm, AppCtx *ctx, DMLabel *adaptLabel)
40568a717dSMatthew G. Knepley {
41568a717dSMatthew G. Knepley   PetscMPIInt    rank;
42568a717dSMatthew G. Knepley 
43568a717dSMatthew G. Knepley   PetscFunctionBegin;
449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank));
459566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Adaptation Label", adaptLabel));
469566063dSJacob Faibussowitsch   if (ctx->refcell[rank] >= 0) PetscCall(DMLabelSetValue(*adaptLabel, ctx->refcell[rank], DM_ADAPT_REFINE));
47568a717dSMatthew G. Knepley   PetscFunctionReturn(0);
48568a717dSMatthew G. Knepley }
49568a717dSMatthew G. Knepley 
50568a717dSMatthew G. Knepley int main(int argc, char **argv)
51568a717dSMatthew G. Knepley {
52568a717dSMatthew G. Knepley   DM             dm, dma;
53568a717dSMatthew G. Knepley   DMLabel        adaptLabel;
54568a717dSMatthew G. Knepley   AppCtx         ctx;
55568a717dSMatthew G. Knepley 
56*327415f7SBarry Smith   PetscFunctionBeginUser;
579566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
589566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx));
599566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm));
609566063dSJacob Faibussowitsch   PetscCall(CreateAdaptLabel(dm, &ctx, &adaptLabel));
619566063dSJacob Faibussowitsch   PetscCall(DMAdaptLabel(dm, adaptLabel, &dma));
629566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) dma, "Adapted Mesh"));
639566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&adaptLabel));
649566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
659566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dma, NULL, "-adapt_dm_view"));
669566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dma));
679566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx.refcell));
689566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
69b122ec5aSJacob Faibussowitsch   return 0;
70568a717dSMatthew G. Knepley }
71568a717dSMatthew G. Knepley 
72568a717dSMatthew G. Knepley /*TEST
73568a717dSMatthew G. Knepley 
74012bc364SMatthew G. Knepley   testset:
75c0517cd5SMatthew G. Knepley     args: -dm_adaptor cellrefiner -dm_plex_transform_type refine_sbr
76012bc364SMatthew G. Knepley 
77568a717dSMatthew G. Knepley     test:
78568a717dSMatthew G. Knepley       suffix: 0
79568a717dSMatthew G. Knepley       requires: triangle
80012bc364SMatthew G. Knepley       args: -dm_view -adapt_dm_view
81568a717dSMatthew G. Knepley 
82568a717dSMatthew G. Knepley     test:
83568a717dSMatthew G. Knepley       suffix: 1
84568a717dSMatthew G. Knepley       requires: triangle
85012bc364SMatthew G. Knepley       args: -dm_coord_space 0 -refcell 2 -dm_view ::ascii_info_detail -adapt_dm_view ::ascii_info_detail
86568a717dSMatthew G. Knepley 
87568a717dSMatthew G. Knepley     test:
88568a717dSMatthew G. Knepley       suffix: 2
89568a717dSMatthew G. Knepley       requires: triangle
90568a717dSMatthew G. Knepley       nsize: 2
91e600fa54SMatthew G. Knepley       args: -refcell 2,-1 -petscpartitioner_type simple -dm_view -adapt_dm_view
92568a717dSMatthew G. Knepley 
93568a717dSMatthew G. Knepley TEST*/
94