1 static char help[] = "Tests for coarsening\n\n"; 2 3 #include <petscdmplex.h> 4 5 typedef struct { 6 PetscBool uninterpolate; /* Uninterpolate the mesh at the end */ 7 } AppCtx; 8 9 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 10 { 11 PetscFunctionBegin; 12 options->uninterpolate = PETSC_FALSE; 13 14 PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX"); 15 PetscCall(PetscOptionsBool("-uninterpolate", "Uninterpolate the mesh at the end", "ex14.c", options->uninterpolate, &options->uninterpolate, NULL)); 16 PetscOptionsEnd(); 17 PetscFunctionReturn(0); 18 } 19 20 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 21 { 22 PetscFunctionBegin; 23 PetscCall(DMCreate(comm, dm)); 24 PetscCall(DMSetType(*dm, DMPLEX)); 25 PetscCall(DMSetFromOptions(*dm)); 26 PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 27 if (user->uninterpolate) { 28 DM udm = NULL; 29 30 PetscCall(DMPlexUninterpolate(*dm, &udm)); 31 PetscCall(DMDestroy(dm)); 32 *dm = udm; 33 PetscCall(DMViewFromOptions(*dm, NULL, "-un_dm_view")); 34 } 35 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 36 PetscFunctionReturn(0); 37 } 38 39 int main(int argc, char **argv) 40 { 41 DM dm; 42 AppCtx user; 43 44 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 45 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user)); 46 PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 47 PetscCall(DMDestroy(&dm)); 48 PetscCall(PetscFinalize()); 49 return 0; 50 } 51 52 /*TEST 53 test: 54 suffix: 0 55 requires: triangle 56 args: -dm_view -dm_refine 1 -dm_coarsen -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_faces 57 58 TEST*/ 59