xref: /petsc/src/dm/impls/plex/tests/ex32.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Tests for periodic mesh output\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscdmplex.h>
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown typedef struct {
6*c4762a1bSJed Brown   DM        dm;
7*c4762a1bSJed Brown   PetscInt  dim;                          /* The topological mesh dimension */
8*c4762a1bSJed Brown   PetscBool cellSimplex;                  /* Use simplices or hexes */
9*c4762a1bSJed Brown   PetscInt  faces[3];                     /* Faces per direction */
10*c4762a1bSJed Brown   PetscBool isPeriodic;                   /* Flag for periodic mesh */
11*c4762a1bSJed Brown   DMBoundaryType periodicity[3];          /* Periodicity per direction */
12*c4762a1bSJed Brown   char      filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */
13*c4762a1bSJed Brown } AppCtx;
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
16*c4762a1bSJed Brown {
17*c4762a1bSJed Brown   PetscInt       n;
18*c4762a1bSJed Brown   PetscErrorCode ierr;
19*c4762a1bSJed Brown 
20*c4762a1bSJed Brown   PetscFunctionBegin;
21*c4762a1bSJed Brown   options->dim            = 2;
22*c4762a1bSJed Brown   options->cellSimplex    = PETSC_TRUE;
23*c4762a1bSJed Brown   options->faces[0]       = 1;
24*c4762a1bSJed Brown   options->faces[1]       = 1;
25*c4762a1bSJed Brown   options->faces[2]       = 1;
26*c4762a1bSJed Brown   options->periodicity[0] = DM_BOUNDARY_NONE;
27*c4762a1bSJed Brown   options->periodicity[1] = DM_BOUNDARY_NONE;
28*c4762a1bSJed Brown   options->periodicity[2] = DM_BOUNDARY_NONE;
29*c4762a1bSJed Brown   options->filename[0]    = '\0';
30*c4762a1bSJed Brown 
31*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex32.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex32.c", options->cellSimplex, &options->cellSimplex, NULL);CHKERRQ(ierr);
34*c4762a1bSJed Brown   n    = 3;
35*c4762a1bSJed Brown   ierr = PetscOptionsIntArray("-faces", "Faces per direction", "ex32.c", options->faces, &n, NULL);CHKERRQ(ierr);
36*c4762a1bSJed Brown   n    = 3;
37*c4762a1bSJed Brown   ierr = PetscOptionsEnumArray("-periodicity", "Periodicity per direction", "ex32.c", DMBoundaryTypes, (PetscEnum *) options->periodicity, &n, &options->isPeriodic);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = PetscOptionsString("-filename", "The mesh file", "ex32.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = PetscOptionsEnd();
40*c4762a1bSJed Brown   PetscFunctionReturn(0);
41*c4762a1bSJed Brown }
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown PetscErrorCode CheckMesh(DM dm, AppCtx *user)
44*c4762a1bSJed Brown {
45*c4762a1bSJed Brown   PetscReal      detJ, J[9], refVol = 1.0;
46*c4762a1bSJed Brown   PetscReal      vol;
47*c4762a1bSJed Brown   PetscInt       dim, depth, d, cStart, cEnd, c;
48*c4762a1bSJed Brown   PetscErrorCode ierr;
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown   PetscFunctionBegin;
51*c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
52*c4762a1bSJed Brown   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
53*c4762a1bSJed Brown   for (d = 0; d < dim; ++d) {
54*c4762a1bSJed Brown     refVol *= 2.0;
55*c4762a1bSJed Brown   }
56*c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
57*c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
58*c4762a1bSJed Brown     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, NULL, J, NULL, &detJ);CHKERRQ(ierr);
59*c4762a1bSJed Brown     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh cell %d is inverted, |J| = %g", c, detJ);
60*c4762a1bSJed Brown     if (depth > 1) {
61*c4762a1bSJed Brown       ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);CHKERRQ(ierr);
62*c4762a1bSJed Brown       if (vol <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mesh cell %d is inverted, vol = %g", c, vol);
63*c4762a1bSJed Brown     }
64*c4762a1bSJed Brown   }
65*c4762a1bSJed Brown   PetscFunctionReturn(0);
66*c4762a1bSJed Brown }
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown PetscErrorCode CompareCones(DM dm, DM idm)
69*c4762a1bSJed Brown {
70*c4762a1bSJed Brown   PetscInt       cStart, cEnd, c, vStart, vEnd, v;
71*c4762a1bSJed Brown   PetscErrorCode ierr;
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown   PetscFunctionBegin;
74*c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
75*c4762a1bSJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
76*c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
77*c4762a1bSJed Brown     const PetscInt *cone;
78*c4762a1bSJed Brown     PetscInt       *points = NULL, numPoints, p, numVertices = 0, coneSize;
79*c4762a1bSJed Brown 
80*c4762a1bSJed Brown     ierr = DMPlexGetCone(dm, c, &cone);CHKERRQ(ierr);
81*c4762a1bSJed Brown     ierr = DMPlexGetConeSize(dm, c, &coneSize);CHKERRQ(ierr);
82*c4762a1bSJed Brown     ierr = DMPlexGetTransitiveClosure(idm, c, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
83*c4762a1bSJed Brown     for (p = 0; p < numPoints*2; p += 2) {
84*c4762a1bSJed Brown       const PetscInt point = points[p];
85*c4762a1bSJed Brown       if ((point >= vStart) && (point < vEnd)) points[numVertices++] = point;
86*c4762a1bSJed Brown     }
87*c4762a1bSJed Brown     if (numVertices != coneSize) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In cell %d, cone size %d != %d vertices in closure", c, coneSize, numVertices);
88*c4762a1bSJed Brown     for (v = 0; v < numVertices; ++v) {
89*c4762a1bSJed Brown       if (cone[v] != points[v]) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "In cell %d, cone point %d is %d != %d vertex in closure", c, v, cone[v], points[v]);
90*c4762a1bSJed Brown     }
91*c4762a1bSJed Brown     ierr = DMPlexRestoreTransitiveClosure(idm, c, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
92*c4762a1bSJed Brown   }
93*c4762a1bSJed Brown   PetscFunctionReturn(0);
94*c4762a1bSJed Brown }
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
97*c4762a1bSJed Brown {
98*c4762a1bSJed Brown   PetscInt       dim         = user->dim;
99*c4762a1bSJed Brown   PetscBool      cellSimplex = user->cellSimplex;
100*c4762a1bSJed Brown   const char    *filename    = user->filename;
101*c4762a1bSJed Brown   size_t         len;
102*c4762a1bSJed Brown   PetscMPIInt    rank;
103*c4762a1bSJed Brown   PetscErrorCode ierr;
104*c4762a1bSJed Brown 
105*c4762a1bSJed Brown   PetscFunctionBegin;
106*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
108*c4762a1bSJed Brown   if (len) {
109*c4762a1bSJed Brown     ierr = DMPlexCreateFromFile(comm, filename, PETSC_FALSE, dm);CHKERRQ(ierr);
110*c4762a1bSJed Brown     ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr);
111*c4762a1bSJed Brown   } else {
112*c4762a1bSJed Brown     PetscReal L[3] = {1.0, 1.0, 1.0};
113*c4762a1bSJed Brown     PetscReal maxCell[3];
114*c4762a1bSJed Brown     PetscInt  d;
115*c4762a1bSJed Brown 
116*c4762a1bSJed Brown     for (d = 0; d < dim; ++d) {maxCell[d] = (1.0/user->faces[d])*1.1;}
117*c4762a1bSJed Brown     ierr = DMPlexCreateBoxMesh(comm, dim, cellSimplex, user->faces, NULL, NULL, user->periodicity, PETSC_TRUE, dm);CHKERRQ(ierr);
118*c4762a1bSJed Brown     ierr = DMSetPeriodicity(*dm, user->isPeriodic, maxCell, L, user->periodicity);CHKERRQ(ierr);
119*c4762a1bSJed Brown   }
120*c4762a1bSJed Brown   {
121*c4762a1bSJed Brown     DM               pdm = NULL;
122*c4762a1bSJed Brown     PetscPartitioner part;
123*c4762a1bSJed Brown 
124*c4762a1bSJed Brown     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
125*c4762a1bSJed Brown     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
126*c4762a1bSJed Brown     /* Distribute mesh over processes */
127*c4762a1bSJed Brown     ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
128*c4762a1bSJed Brown     if (pdm) {
129*c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
130*c4762a1bSJed Brown       *dm  = pdm;
131*c4762a1bSJed Brown     }
132*c4762a1bSJed Brown   }
133*c4762a1bSJed Brown   ierr = DMLocalizeCoordinates(*dm);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Mesh");CHKERRQ(ierr);
135*c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
136*c4762a1bSJed Brown   user->dm = *dm;
137*c4762a1bSJed Brown   PetscFunctionReturn(0);
138*c4762a1bSJed Brown }
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown int main(int argc, char **argv)
141*c4762a1bSJed Brown {
142*c4762a1bSJed Brown   AppCtx         user;                 /* user-defined work context */
143*c4762a1bSJed Brown   PetscErrorCode ierr;
144*c4762a1bSJed Brown 
145*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
146*c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &user.dm);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = CheckMesh(user.dm, &user);CHKERRQ(ierr);
149*c4762a1bSJed Brown   ierr = DMDestroy(&user.dm);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = PetscFinalize();
151*c4762a1bSJed Brown   return ierr;
152*c4762a1bSJed Brown }
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown /*TEST
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown   test:
157*c4762a1bSJed Brown     suffix: 0
158*c4762a1bSJed Brown     args: -dim 2 -cell_simplex 0 -faces 3,1,0 -periodicity periodic,none,none -dm_view ::ascii_info_detail
159*c4762a1bSJed Brown   test:
160*c4762a1bSJed Brown     suffix: 1
161*c4762a1bSJed Brown     nsize: 2
162*c4762a1bSJed Brown     args: -dim 2 -cell_simplex 0 -faces 3,1,0 -periodicity periodic,none,none -petscpartitioner_type simple -dm_view ::ascii_info_detail
163*c4762a1bSJed Brown   test:
164*c4762a1bSJed Brown     suffix: 2
165*c4762a1bSJed Brown     nsize: 2
166*c4762a1bSJed Brown     args: -dim 2 -cell_simplex 0 -faces 6,2,0 -periodicity periodic,none,none -petscpartitioner_type simple -dm_view ::ascii_info_detail
167*c4762a1bSJed Brown   test:
168*c4762a1bSJed Brown     suffix: 3
169*c4762a1bSJed Brown     nsize: 4
170*c4762a1bSJed Brown     args: -dim 2 -cell_simplex 0 -faces 6,2,0 -periodicity periodic,none,none -petscpartitioner_type simple -dm_view ::ascii_info_detail
171*c4762a1bSJed Brown   test:
172*c4762a1bSJed Brown     suffix: 4
173*c4762a1bSJed Brown     nsize: 2
174*c4762a1bSJed Brown     args: -dim 2 -cell_simplex 0 -faces 3,1,0 -periodicity periodic,none,none -dm_plex_periodic_cut -petscpartitioner_type simple -dm_view ::ascii_info_detail
175*c4762a1bSJed Brown   test:
176*c4762a1bSJed Brown     suffix: 5
177*c4762a1bSJed Brown     nsize: 2
178*c4762a1bSJed Brown     args: -dim 2 -cell_simplex 0 -faces 6,2,0 -periodicity periodic,none,none -dm_plex_periodic_cut -petscpartitioner_type simple -dm_view ::ascii_info_detail
179*c4762a1bSJed Brown   test:
180*c4762a1bSJed Brown     suffix: 6
181*c4762a1bSJed Brown     nsize: 4
182*c4762a1bSJed Brown     args: -dim 2 -cell_simplex 0 -faces 6,2,0 -periodicity periodic,none,none -dm_plex_periodic_cut -petscpartitioner_type simple -dm_view ::ascii_info_detail
183*c4762a1bSJed Brown 
184*c4762a1bSJed Brown TEST*/
185