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