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 { 8 DM dm, dmAdapt; 9 DMLabel adaptLabel; 10 PetscInt cStart, cEnd; 11 12 PetscFunctionBeginUser; 13 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 14 PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 15 PetscCall(DMSetType(dm, DMPLEX)); 16 PetscCall(DMSetFromOptions(dm)); 17 PetscCall(PetscObjectSetName((PetscObject) dm, "Pre Adaptation Mesh")); 18 PetscCall(DMViewFromOptions(dm, NULL, "-pre_adapt_dm_view")); 19 20 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 21 PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 22 PetscCall(DMLabelSetDefaultValue(adaptLabel, DM_ADAPT_COARSEN)); 23 if (cEnd > cStart) PetscCall(DMLabelSetValue(adaptLabel, cStart, DM_ADAPT_REFINE)); 24 PetscCall(DMAdaptLabel(dm, adaptLabel, &dmAdapt)); 25 PetscCall(PetscObjectSetName((PetscObject) dmAdapt, "Post Adaptation Mesh")); 26 PetscCall(DMViewFromOptions(dmAdapt, NULL, "-post_adapt_dm_view")); 27 PetscCall(DMDestroy(&dmAdapt)); 28 PetscCall(DMLabelDestroy(&adaptLabel)); 29 PetscCall(DMDestroy(&dm)); 30 PetscCall(PetscFinalize()); 31 return 0; 32 } 33 34 /*TEST 35 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 40 test: 41 suffix: 3d_tetgen 42 requires: tetgen 43 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 44 test: 45 suffix: 3d_ctetgen 46 requires: ctetgen !complex !single 47 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 48 49 TEST*/ 50