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