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