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