xref: /petsc/src/dm/impls/plex/tests/ex20.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown const char help[] = "Test DMPlex implementation of DMAdaptLabel().\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdm.h>
4c4762a1bSJed Brown #include <petscdmplex.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc, char **argv)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   DM             dm, dmAdapt;
9c4762a1bSJed Brown   DMLabel        adaptLabel;
1030602db0SMatthew G. Knepley   PetscInt       cStart, cEnd;
11c4762a1bSJed Brown 
12*327415f7SBarry Smith   PetscFunctionBeginUser;
139566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
149566063dSJacob Faibussowitsch   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
159566063dSJacob Faibussowitsch   PetscCall(DMSetType(dm, DMPLEX));
169566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
179566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) dm, "Pre Adaptation Mesh"));
189566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm, NULL, "-pre_adapt_dm_view"));
1930602db0SMatthew G. Knepley 
209566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
219566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
229566063dSJacob Faibussowitsch   PetscCall(DMLabelSetDefaultValue(adaptLabel, DM_ADAPT_COARSEN));
239566063dSJacob Faibussowitsch   if (cEnd > cStart) PetscCall(DMLabelSetValue(adaptLabel, cStart, DM_ADAPT_REFINE));
249566063dSJacob Faibussowitsch   PetscCall(DMAdaptLabel(dm, adaptLabel, &dmAdapt));
259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) dmAdapt, "Post Adaptation Mesh"));
269566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dmAdapt, NULL, "-post_adapt_dm_view"));
279566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dmAdapt));
289566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&adaptLabel));
299566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
309566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
31b122ec5aSJacob Faibussowitsch   return 0;
32c4762a1bSJed Brown }
33c4762a1bSJed Brown 
34c4762a1bSJed Brown /*TEST
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   test:
37c4762a1bSJed Brown     suffix: 2d
38c4762a1bSJed Brown     requires: triangle !single
391c35d092SJoe Wallwork     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
40c4762a1bSJed Brown   test:
41c4762a1bSJed Brown     suffix: 3d_tetgen
42b5a892a1SMatthew G. Knepley     requires: tetgen
431c35d092SJoe Wallwork     args: -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_coord_space 0 -pre_adapt_dm_view ascii::ascii_info -post_adapt_dm_view ascii::ascii_info
44c4762a1bSJed Brown   test:
45c4762a1bSJed Brown     suffix: 3d_ctetgen
46c4762a1bSJed Brown     requires: ctetgen !complex !single
471c35d092SJoe Wallwork     args: -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_coord_space 0 -pre_adapt_dm_view ascii::ascii_info -post_adapt_dm_view ascii::ascii_info
48c4762a1bSJed Brown 
49c4762a1bSJed Brown TEST*/
50