xref: /petsc/src/dm/impls/plex/tests/ex14.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008) !
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   PetscFunctionBeginUser;
45   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
46   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
47   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
48   PetscCall(DMDestroy(&dm));
49   PetscCall(PetscFinalize());
50   return 0;
51 }
52 
53 /*TEST
54   test:
55     suffix: 0
56     requires: triangle
57     args: -dm_view -dm_refine 1 -dm_coarsen -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_faces
58 
59 TEST*/
60