xref: /petsc/src/dm/impls/plex/tests/ex41.c (revision 568a717df64ec89e8f259c8e7378de26d0a21de1)
1*568a717dSMatthew G. Knepley static const char help[] = "Tests for adaptive refinement";
2*568a717dSMatthew G. Knepley 
3*568a717dSMatthew G. Knepley #include <petscdmplex.h>
4*568a717dSMatthew G. Knepley 
5*568a717dSMatthew G. Knepley typedef struct {
6*568a717dSMatthew G. Knepley   PetscBool metric;  /* Flag to use metric adaptation, instead of tagging */
7*568a717dSMatthew G. Knepley   PetscInt *refcell; /* A cell to be refined on each process */
8*568a717dSMatthew G. Knepley } AppCtx;
9*568a717dSMatthew G. Knepley 
10*568a717dSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
11*568a717dSMatthew G. Knepley {
12*568a717dSMatthew G. Knepley   PetscMPIInt    size;
13*568a717dSMatthew G. Knepley   PetscInt       n;
14*568a717dSMatthew G. Knepley   PetscErrorCode ierr;
15*568a717dSMatthew G. Knepley 
16*568a717dSMatthew G. Knepley   PetscFunctionBeginUser;
17*568a717dSMatthew G. Knepley   options->metric  = PETSC_FALSE;
18*568a717dSMatthew G. Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
19*568a717dSMatthew G. Knepley   ierr = PetscCalloc1(size, &options->refcell);CHKERRQ(ierr);
20*568a717dSMatthew G. Knepley   n    = size;
21*568a717dSMatthew G. Knepley 
22*568a717dSMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "Parallel Mesh Adaptation Options", "DMPLEX");CHKERRQ(ierr);
23*568a717dSMatthew G. Knepley   ierr = PetscOptionsBool("-metric", "Flag for metric refinement", "ex41.c", options->metric, &options->metric, NULL);CHKERRQ(ierr);
24*568a717dSMatthew G. Knepley   ierr = PetscOptionsIntArray("-refcell", "The cell to be refined", "ex41.c", options->refcell, &n, NULL);CHKERRQ(ierr);
25*568a717dSMatthew G. Knepley   if (n && n != size) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Only gave %D cells to refine, must give one for all %D processes", n, size);
26*568a717dSMatthew G. Knepley   ierr = PetscOptionsEnd();
27*568a717dSMatthew G. Knepley   PetscFunctionReturn(0);
28*568a717dSMatthew G. Knepley }
29*568a717dSMatthew G. Knepley 
30*568a717dSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
31*568a717dSMatthew G. Knepley {
32*568a717dSMatthew G. Knepley   PetscErrorCode ierr;
33*568a717dSMatthew G. Knepley 
34*568a717dSMatthew G. Knepley   PetscFunctionBegin;
35*568a717dSMatthew G. Knepley   ierr = DMPlexCreateBoxMesh(comm, 2, PETSC_TRUE, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
36*568a717dSMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
37*568a717dSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
38*568a717dSMatthew G. Knepley   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
39*568a717dSMatthew G. Knepley   PetscFunctionReturn(0);
40*568a717dSMatthew G. Knepley }
41*568a717dSMatthew G. Knepley 
42*568a717dSMatthew G. Knepley static PetscErrorCode CreateAdaptLabel(DM dm, AppCtx *ctx, DMLabel *adaptLabel)
43*568a717dSMatthew G. Knepley {
44*568a717dSMatthew G. Knepley   PetscMPIInt    rank;
45*568a717dSMatthew G. Knepley   PetscErrorCode ierr;
46*568a717dSMatthew G. Knepley 
47*568a717dSMatthew G. Knepley   PetscFunctionBegin;
48*568a717dSMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
49*568a717dSMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, "Adaptation Label", adaptLabel);CHKERRQ(ierr);
50*568a717dSMatthew G. Knepley   if (ctx->refcell[rank] >= 0) {ierr = DMLabelSetValue(*adaptLabel, ctx->refcell[rank], DM_ADAPT_REFINE);CHKERRQ(ierr);}
51*568a717dSMatthew G. Knepley   PetscFunctionReturn(0);
52*568a717dSMatthew G. Knepley }
53*568a717dSMatthew G. Knepley 
54*568a717dSMatthew G. Knepley int main(int argc, char **argv)
55*568a717dSMatthew G. Knepley {
56*568a717dSMatthew G. Knepley   DM             dm, dma;
57*568a717dSMatthew G. Knepley   DMLabel        adaptLabel;
58*568a717dSMatthew G. Knepley   AppCtx         ctx;
59*568a717dSMatthew G. Knepley   PetscErrorCode ierr;
60*568a717dSMatthew G. Knepley 
61*568a717dSMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr;
62*568a717dSMatthew G. Knepley   ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr);
63*568a717dSMatthew G. Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &ctx, &dm);CHKERRQ(ierr);
64*568a717dSMatthew G. Knepley   ierr = CreateAdaptLabel(dm, &ctx, &adaptLabel);CHKERRQ(ierr);
65*568a717dSMatthew G. Knepley   ierr = DMAdaptLabel(dm, adaptLabel, &dma);CHKERRQ(ierr);
66*568a717dSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) dma, "Adapted Mesh");CHKERRQ(ierr);
67*568a717dSMatthew G. Knepley   ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr);
68*568a717dSMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
69*568a717dSMatthew G. Knepley   ierr = DMViewFromOptions(dma, NULL, "-adapt_dm_view");CHKERRQ(ierr);
70*568a717dSMatthew G. Knepley   ierr = DMDestroy(&dma);CHKERRQ(ierr);
71*568a717dSMatthew G. Knepley   ierr = PetscFree(ctx.refcell);CHKERRQ(ierr);
72*568a717dSMatthew G. Knepley   ierr = PetscFinalize();
73*568a717dSMatthew G. Knepley   return ierr;
74*568a717dSMatthew G. Knepley }
75*568a717dSMatthew G. Knepley 
76*568a717dSMatthew G. Knepley /*TEST
77*568a717dSMatthew G. Knepley 
78*568a717dSMatthew G. Knepley   test:
79*568a717dSMatthew G. Knepley     suffix: 0
80*568a717dSMatthew G. Knepley     requires: triangle
81*568a717dSMatthew G. Knepley     args: -dm_plex_adaptor cellrefiner -dm_plex_cell_refiner sbr -dm_view -adapt_dm_view
82*568a717dSMatthew G. Knepley 
83*568a717dSMatthew G. Knepley   test:
84*568a717dSMatthew G. Knepley     suffix: 1
85*568a717dSMatthew G. Knepley     requires: triangle
86*568a717dSMatthew G. Knepley     args: -refcell 2 -dm_plex_adaptor cellrefiner -dm_plex_cell_refiner sbr -dm_view ::ascii_info_detail -adapt_dm_view ::ascii_info_detail
87*568a717dSMatthew G. Knepley 
88*568a717dSMatthew G. Knepley   test:
89*568a717dSMatthew G. Knepley     suffix: 2
90*568a717dSMatthew G. Knepley     requires: triangle
91*568a717dSMatthew G. Knepley     nsize: 2
92*568a717dSMatthew G. Knepley     args: -refcell 2,-1 -dm_distribute -petscpartitioner_type simple -dm_plex_adaptor cellrefiner -dm_plex_cell_refiner sbr -dm_view -adapt_dm_view
93*568a717dSMatthew G. Knepley 
94*568a717dSMatthew G. Knepley TEST*/
95