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 PetscErrorCode ierr; 12 13 PetscFunctionBegin; 14 options->uninterpolate = PETSC_FALSE; 15 16 ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr); 17 CHKERRQ(PetscOptionsBool("-uninterpolate", "Uninterpolate the mesh at the end", "ex14.c", options->uninterpolate, &options->uninterpolate, NULL)); 18 ierr = PetscOptionsEnd();CHKERRQ(ierr); 19 PetscFunctionReturn(0); 20 } 21 22 PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) 23 { 24 PetscFunctionBegin; 25 CHKERRQ(DMCreate(comm, dm)); 26 CHKERRQ(DMSetType(*dm, DMPLEX)); 27 CHKERRQ(DMSetFromOptions(*dm)); 28 CHKERRQ(DMViewFromOptions(*dm, NULL, "-orig_dm_view")); 29 if (user->uninterpolate) { 30 DM udm = NULL; 31 32 CHKERRQ(DMPlexUninterpolate(*dm, &udm)); 33 CHKERRQ(DMDestroy(dm)); 34 *dm = udm; 35 CHKERRQ(DMViewFromOptions(*dm, NULL, "-un_dm_view")); 36 } 37 CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view")); 38 PetscFunctionReturn(0); 39 } 40 41 int main(int argc, char **argv) 42 { 43 DM dm; 44 AppCtx user; 45 PetscErrorCode ierr; 46 47 ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr; 48 CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user)); 49 CHKERRQ(CreateMesh(PETSC_COMM_WORLD, &user, &dm)); 50 CHKERRQ(DMDestroy(&dm)); 51 ierr = PetscFinalize(); 52 return ierr; 53 } 54 55 /*TEST 56 test: 57 suffix: 0 58 requires: triangle 59 args: -dm_view -dm_refine 1 -dm_coarsen -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_faces 60 61 TEST*/ 62