xref: /petsc/src/dm/impls/plex/tests/ex20.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 const char help[] = "Test DMPlex implementation of DMAdaptLabel().\n\n";
2 
3 #include <petscdm.h>
4 #include <petscdmplex.h>
5 
6 int main(int argc, char **argv) {
7   DM       dm, dmAdapt;
8   DMLabel  adaptLabel;
9   PetscInt cStart, cEnd;
10 
11   PetscFunctionBeginUser;
12   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
13   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
14   PetscCall(DMSetType(dm, DMPLEX));
15   PetscCall(DMSetFromOptions(dm));
16   PetscCall(PetscObjectSetName((PetscObject)dm, "Pre Adaptation Mesh"));
17   PetscCall(DMViewFromOptions(dm, NULL, "-pre_adapt_dm_view"));
18 
19   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
20   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
21   PetscCall(DMLabelSetDefaultValue(adaptLabel, DM_ADAPT_COARSEN));
22   if (cEnd > cStart) PetscCall(DMLabelSetValue(adaptLabel, cStart, DM_ADAPT_REFINE));
23   PetscCall(DMAdaptLabel(dm, adaptLabel, &dmAdapt));
24   PetscCall(PetscObjectSetName((PetscObject)dmAdapt, "Post Adaptation Mesh"));
25   PetscCall(DMViewFromOptions(dmAdapt, NULL, "-post_adapt_dm_view"));
26   PetscCall(DMDestroy(&dmAdapt));
27   PetscCall(DMLabelDestroy(&adaptLabel));
28   PetscCall(DMDestroy(&dm));
29   PetscCall(PetscFinalize());
30   return 0;
31 }
32 
33 /*TEST
34 
35   test:
36     suffix: 2d
37     requires: triangle !single
38     args: -dm_plex_box_faces 3,3 -dm_coord_space 0 -pre_adapt_dm_view ascii::ascii_info -post_adapt_dm_view ascii::ascii_info
39 
40   # We eliminate the lines with "marker" because different compiler flags make the meshes produce different surface meshes
41   testset:
42     args: -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_coord_space 0 \
43           -pre_adapt_dm_view ascii::ascii_info -post_adapt_dm_view ascii::ascii_info
44     filter: grep -v "marker"
45 
46     test:
47       suffix: 3d_tetgen
48       requires: tetgen
49 
50     test:
51       suffix: 3d_ctetgen
52       requires: ctetgen !complex !single
53 
54 TEST*/
55