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