xref: /petsc/src/dm/impls/plex/tests/ex16.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Tests for creation of submeshes\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscdmplex.h>
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown typedef struct {
6*c4762a1bSJed Brown   PetscInt  debug;       /* The debugging level */
7*c4762a1bSJed Brown   PetscInt  dim;         /* The topological mesh dimension */
8*c4762a1bSJed Brown   PetscBool cellSimplex; /* Use simplices or hexes */
9*c4762a1bSJed Brown } AppCtx;
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12*c4762a1bSJed Brown {
13*c4762a1bSJed Brown   PetscErrorCode ierr;
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown   PetscFunctionBegin;
16*c4762a1bSJed Brown   options->debug       = 0;
17*c4762a1bSJed Brown   options->dim         = 2;
18*c4762a1bSJed Brown   options->cellSimplex = PETSC_TRUE;
19*c4762a1bSJed Brown 
20*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
21*c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-debug", "The debugging level", "ex16.c", options->debug, &options->debug, NULL,0);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex16.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
23*c4762a1bSJed Brown   ierr = PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex16.c", options->cellSimplex, &options->cellSimplex, NULL);CHKERRQ(ierr);
24*c4762a1bSJed Brown   ierr = PetscOptionsEnd();
25*c4762a1bSJed Brown   PetscFunctionReturn(0);
26*c4762a1bSJed Brown }
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
29*c4762a1bSJed Brown {
30*c4762a1bSJed Brown   PetscInt       dim         = user->dim;
31*c4762a1bSJed Brown   PetscBool      cellSimplex = user->cellSimplex;
32*c4762a1bSJed Brown   PetscErrorCode ierr;
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown   PetscFunctionBegin;
35*c4762a1bSJed Brown   ierr = DMPlexCreateBoxMesh(comm, dim, cellSimplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
36*c4762a1bSJed Brown   {
37*c4762a1bSJed Brown     DM pdm = NULL;
38*c4762a1bSJed Brown 
39*c4762a1bSJed Brown     /* Distribute mesh over processes */
40*c4762a1bSJed Brown     ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
41*c4762a1bSJed Brown     if (pdm) {
42*c4762a1bSJed Brown       ierr = DMViewFromOptions(pdm, NULL, "-dm_view");CHKERRQ(ierr);
43*c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
44*c4762a1bSJed Brown       *dm  = pdm;
45*c4762a1bSJed Brown     }
46*c4762a1bSJed Brown   }
47*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
50*c4762a1bSJed Brown   PetscFunctionReturn(0);
51*c4762a1bSJed Brown }
52*c4762a1bSJed Brown 
53*c4762a1bSJed Brown PetscErrorCode CreateSubmesh(DM dm, PetscBool start, DM *subdm)
54*c4762a1bSJed Brown {
55*c4762a1bSJed Brown   DMLabel        label, map;
56*c4762a1bSJed Brown   PetscInt       cStart, cEnd, cStartSub, cEndSub, c;
57*c4762a1bSJed Brown   PetscErrorCode ierr;
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown   PetscFunctionBegin;
60*c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
61*c4762a1bSJed Brown   ierr = DMCreateLabel(dm, "cells");CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = DMGetLabel(dm, "cells", &label);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = DMLabelClearStratum(label, 1);CHKERRQ(ierr);
64*c4762a1bSJed Brown   if (start) {cStartSub = cStart; cEndSub = cEnd/2;}
65*c4762a1bSJed Brown   else       {cStartSub = cEnd/2; cEndSub = cEnd;}
66*c4762a1bSJed Brown   for (c = cStartSub; c < cEndSub; ++c) {ierr = DMLabelSetValue(label, c, 1);CHKERRQ(ierr);}
67*c4762a1bSJed Brown   ierr = DMPlexFilter(dm, label, 1, subdm);CHKERRQ(ierr);
68*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *subdm, "Submesh");CHKERRQ(ierr);
69*c4762a1bSJed Brown   ierr = DMViewFromOptions(*subdm, NULL, "-dm_view");CHKERRQ(ierr);
70*c4762a1bSJed Brown   ierr = DMPlexGetSubpointMap(*subdm, &map);CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = DMLabelView(map, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
72*c4762a1bSJed Brown   PetscFunctionReturn(0);
73*c4762a1bSJed Brown }
74*c4762a1bSJed Brown 
75*c4762a1bSJed Brown int main(int argc, char **argv)
76*c4762a1bSJed Brown {
77*c4762a1bSJed Brown   DM             dm, subdm;
78*c4762a1bSJed Brown   AppCtx         user;
79*c4762a1bSJed Brown   PetscErrorCode ierr;
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
82*c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = CreateSubmesh(dm, PETSC_TRUE, &subdm);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = DMSetFromOptions(subdm);CHKERRQ(ierr);
86*c4762a1bSJed Brown   ierr = DMDestroy(&subdm);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = CreateSubmesh(dm, PETSC_FALSE, &subdm);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = DMSetFromOptions(subdm);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = DMDestroy(&subdm);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = PetscFinalize();
92*c4762a1bSJed Brown   return ierr;
93*c4762a1bSJed Brown }
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown /*TEST
96*c4762a1bSJed Brown 
97*c4762a1bSJed Brown   test:
98*c4762a1bSJed Brown     suffix: 0
99*c4762a1bSJed Brown     requires: triangle
100*c4762a1bSJed Brown     args: -dm_view ascii::ascii_info_detail -dm_plex_check_symmetry -dm_plex_check_skeleton -dm_plex_check_faces
101*c4762a1bSJed Brown 
102*c4762a1bSJed Brown TEST*/
103