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