xref: /petsc/src/dm/impls/plex/tests/ex20.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(PETSC_COMM_WORLD, &dm));
145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(dm, DMPLEX));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(dm));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) dm, "Pre Adaptation Mesh"));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm, NULL, "-pre_adapt_dm_view"));
1830602db0SMatthew G. Knepley 
195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelSetDefaultValue(adaptLabel, DM_ADAPT_COARSEN));
225f80ce2aSJacob Faibussowitsch   if (cEnd > cStart) CHKERRQ(DMLabelSetValue(adaptLabel, cStart, DM_ADAPT_REFINE));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(DMAdaptLabel(dm, adaptLabel, &dmAdapt));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) dmAdapt, "Post Adaptation Mesh"));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dmAdapt, NULL, "-post_adapt_dm_view"));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dmAdapt));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLabelDestroy(&adaptLabel));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
29*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
30*b122ec5aSJacob Faibussowitsch   return 0;
31c4762a1bSJed Brown }
32c4762a1bSJed Brown 
33c4762a1bSJed Brown /*TEST
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   test:
36c4762a1bSJed Brown     suffix: 2d
37c4762a1bSJed Brown     requires: triangle !single
381c35d092SJoe 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
39c4762a1bSJed Brown   test:
40c4762a1bSJed Brown     suffix: 3d_tetgen
41b5a892a1SMatthew G. Knepley     requires: tetgen
421c35d092SJoe 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
43c4762a1bSJed Brown   test:
44c4762a1bSJed Brown     suffix: 3d_ctetgen
45c4762a1bSJed Brown     requires: ctetgen !complex !single
461c35d092SJoe 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
47c4762a1bSJed Brown 
48c4762a1bSJed Brown TEST*/
49