DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
Hex8.hpp
Go to the documentation of this file.
1#pragma once
2// Auto-generated by tools/gen_shape_functions -- DO NOT EDIT
3// Element: Hex8
4// Regenerate: /usr/bin/python3 -m tools.gen_shape_functions.generate
5
6#include "DNDS/Defines.hpp"
7#include "Geom/Geometric.hpp"
8#include "Geom/ElemEnum.hpp"
10
11namespace DNDS::Geom::Elem
12{
13
14 // Forward declaration (primary template is in ElementTraitsBase.hpp)
15 template <ElemType> struct ShapeFuncImpl;
16
17 // <GEN_SHAPE_FUNCS_BEGIN>
18 template <>
20 {
21 template <class TPoint, class TArray>
22 DNDS_DEVICE_CALLABLE static inline void Diff0(const TPoint &p, TArray &&v)
23 {
24 t_real xi = p[0];
25 t_real et = p[1];
26 t_real zt = p[2];
27 const t_real _t0 = xi - 1;
28 const t_real _t1 = et - 1;
29 const t_real _t2 = ((0.125))*zt + (-0.125);
30 const t_real _t3 = _t1*_t2;
31 const t_real _t4 = xi + 1;
32 const t_real _t5 = et + 1;
33 const t_real _t6 = _t2*_t5;
34 const t_real _t7 = ((0.125))*zt + (0.125);
35 const t_real _t8 = _t1*_t7;
36 const t_real _t9 = _t5*_t7;
37 v(0, 0) = -_t0*_t3;
38 v(0, 1) = _t3*_t4;
39 v(0, 2) = -_t4*_t6;
40 v(0, 3) = _t0*_t6;
41 v(0, 4) = _t0*_t8;
42 v(0, 5) = -_t4*_t8;
43 v(0, 6) = _t4*_t9;
44 v(0, 7) = -_t0*_t9;
45 }
46
47 template <class TPoint, class TArray>
48 DNDS_DEVICE_CALLABLE static inline void Diff1(const TPoint &p, TArray &&v)
49 {
50 t_real xi = p[0];
51 t_real et = p[1];
52 t_real zt = p[2];
53 const t_real _t0 = et - 1;
54 const t_real _t1 = ((0.125))*zt + (-0.125);
55 const t_real _t2 = _t0*_t1;
56 const t_real _t3 = et + 1;
57 const t_real _t4 = _t1*_t3;
58 const t_real _t5 = ((0.125))*zt + (0.125);
59 const t_real _t6 = _t0*_t5;
60 const t_real _t7 = _t3*_t5;
61 const t_real _t8 = xi - 1;
62 const t_real _t9 = _t1*_t8;
63 const t_real _t10 = xi + 1;
64 const t_real _t11 = _t1*_t10;
65 const t_real _t12 = _t5*_t8;
66 const t_real _t13 = _t10*_t5;
67 const t_real _t14 = ((0.125))*_t0;
68 const t_real _t15 = _t14*_t8;
69 const t_real _t16 = _t10*_t14;
70 const t_real _t17 = ((0.125))*_t3;
71 const t_real _t18 = _t10*_t17;
72 const t_real _t19 = _t17*_t8;
73 v(0, 0) = -_t2;
74 v(0, 1) = _t2;
75 v(0, 2) = -_t4;
76 v(0, 3) = _t4;
77 v(0, 4) = _t6;
78 v(0, 5) = -_t6;
79 v(0, 6) = _t7;
80 v(0, 7) = -_t7;
81 v(1, 0) = -_t9;
82 v(1, 1) = _t11;
83 v(1, 2) = -_t11;
84 v(1, 3) = _t9;
85 v(1, 4) = _t12;
86 v(1, 5) = -_t13;
87 v(1, 6) = _t13;
88 v(1, 7) = -_t12;
89 v(2, 0) = -_t15;
90 v(2, 1) = _t16;
91 v(2, 2) = -_t18;
92 v(2, 3) = _t19;
93 v(2, 4) = _t15;
94 v(2, 5) = -_t16;
95 v(2, 6) = _t18;
96 v(2, 7) = -_t19;
97 }
98
99 template <class TPoint, class TArray>
100 DNDS_DEVICE_CALLABLE static inline void Diff2(const TPoint &p, TArray &&v)
101 {
102 t_real xi = p[0];
103 t_real et = p[1];
104 t_real zt = p[2];
105 const t_real _t0 = ((0.125))*zt;
106 const t_real _t1 = _t0 + (-0.125);
107 const t_real _t2 = -_t1;
108 const t_real _t3 = _t0 + (0.125);
109 const t_real _t4 = -_t3;
110 const t_real _t5 = ((0.125))*xi;
111 const t_real _t6 = _t5 + (-0.125);
112 const t_real _t7 = -_t6;
113 const t_real _t8 = _t5 + (0.125);
114 const t_real _t9 = -_t8;
115 const t_real _t10 = ((0.125))*et;
116 const t_real _t11 = _t10 + (-0.125);
117 const t_real _t12 = -_t11;
118 const t_real _t13 = _t10 + (0.125);
119 const t_real _t14 = -_t13;
120 v(3, 0) = _t2;
121 v(3, 1) = _t1;
122 v(3, 2) = _t2;
123 v(3, 3) = _t1;
124 v(3, 4) = _t3;
125 v(3, 5) = _t4;
126 v(3, 6) = _t3;
127 v(3, 7) = _t4;
128 v(4, 0) = _t7;
129 v(4, 1) = _t8;
130 v(4, 2) = _t9;
131 v(4, 3) = _t6;
132 v(4, 4) = _t6;
133 v(4, 5) = _t9;
134 v(4, 6) = _t8;
135 v(4, 7) = _t7;
136 v(5, 0) = _t12;
137 v(5, 1) = _t11;
138 v(5, 2) = _t14;
139 v(5, 3) = _t13;
140 v(5, 4) = _t11;
141 v(5, 5) = _t12;
142 v(5, 6) = _t13;
143 v(5, 7) = _t14;
144 }
145
146 template <class TPoint, class TArray>
147 DNDS_DEVICE_CALLABLE static inline void Diff3(const TPoint &p, TArray &&v)
148 {
149 t_real xi = p[0];
150 t_real et = p[1];
151 t_real zt = p[2];
152 v(9, 0) = (-0.125);
153 v(9, 1) = (0.125);
154 v(9, 2) = (-0.125);
155 v(9, 3) = (0.125);
156 v(9, 4) = (0.125);
157 v(9, 5) = (-0.125);
158 v(9, 6) = (0.125);
159 v(9, 7) = (-0.125);
160 }
161 };
162 // <GEN_SHAPE_FUNCS_END>
163
164
165 /**
166 * @brief Element traits for 8-node trilinear hexahedron (Hex8)
167 *
168 * Hex8 is the simplest 3D hexahedral element, commonly used for:
169 * - Structured hexahedral meshes
170 * - Boundary-fitted grids
171 * - Isoparametric mapping of complex geometries
172 *
173 * Geometry:
174 * - Reference hex: xi, eta, zeta in [-1, 1]
175 * - 8 nodes at the corners of the reference cube
176 * - Parametric space volume = 8.0
177 *
178 * Faces:
179 * - 6 faces, each is a Quad4 element
180 */
181 template <>
183 {
184 // ============================================================
185 // Core Element Identification
186 // ============================================================
187
188 static constexpr ElemType elemType = Hex8;
189 static constexpr int dim = 3;
190 static constexpr int order = 1;
191 static constexpr int numVertices = 8;
192 static constexpr int numNodes = 8;
193 static constexpr int numFaces = 6;
194 static constexpr ParamSpace paramSpace = HexSpace;
195 static constexpr t_real paramSpaceVol = 8.0;
196
197 // ============================================================
198 // Geometry Definition
199 // ============================================================
200
201 /**
202 * @brief Standard coordinates of nodes in parametric space
203 *
204 * Reference cube [-1,1]^3 with nodes at corners.
205 * Standard ordering: bottom face (z=-1) counter-clockwise,
206 * then top face (z=+1) counter-clockwise.
207 *
208 * Node 0-3: bottom face (z = -1)
209 * Node 4-7: top face (z = +1)
210 */
211 static constexpr std::array<t_real, 3 * 8> standardCoords = {
212 -1, -1, -1, // Node 0: bottom-back-left
213 1, -1, -1, // Node 1: bottom-back-right
214 1, 1, -1, // Node 2: bottom-front-right
215 -1, 1, -1, // Node 3: bottom-front-left
216 -1, -1, 1, // Node 4: top-back-left
217 1, -1, 1, // Node 5: top-back-right
218 1, 1, 1, // Node 6: top-front-right
219 -1, 1, 1}; // Node 7: top-front-left
220
221 // ============================================================
222 // Face/Edge Definitions
223 // ============================================================
224
225 /**
226 * @brief Get the element type of a face
227 * @return Quad4 (all faces of hex are bilinear quads)
228 */
229 static constexpr ElemType GetFaceType(t_index /*iFace*/) { return Quad4; }
230
231 /**
232 * @brief Node indices for each face (bilinear quad)
233 *
234 * Face 0: bottom face (z=-1, nodes 0-3-2-1) - note reversed for outward normal
235 * Face 1: back face (y=-1, nodes 0-1-5-4)
236 * Face 2: right face (x=+1, nodes 1-2-6-5)
237 * Face 3: front face (y=+1, nodes 2-3-7-6)
238 * Face 4: left face (x=-1, nodes 3-0-4-7)
239 * Face 5: top face (z=+1, nodes 4-5-6-7)
240 */
241 static constexpr std::array<std::array<t_index, 10>, 6> faceNodes = {{
242 {0, 3, 2, 1}, // Face 0: bottom (z=-1)
243 {0, 1, 5, 4}, // Face 1: back (y=-1)
244 {1, 2, 6, 5}, // Face 2: right (x=+1)
245 {2, 3, 7, 6}, // Face 3: front (y=+1)
246 {0, 4, 7, 3}, // Face 4: left (x=-1)
247 {4, 5, 6, 7}}}; // Face 5: top (z=+1)
248
249 // ============================================================
250 // Order Elevation (P-Refinement)
251 // ============================================================
252
253 /**
254 * @brief Element type after order elevation (O1 -> O2)
255 * Hex8 elevates to Hex27 (27-node triquadratic hex)
256 */
257 static constexpr ElemType elevatedType = Hex27;
258
259 /// @brief Number of additional nodes created during elevation
260 static constexpr int numElevNodes = 19;
261
262 /**
263 * @brief Elevation spans define new node connections
264 *
265 * Elevation creates 19 new nodes:
266 * Spans 0-11: edge midpoints (12 edges)
267 * Spans 12-17: face centers (6 faces)
268 * Span 18: body center (all 8 vertices)
269 */
270 static constexpr std::array<tElevSpan, 19> elevSpans = {{
271 // Bottom face edges
272 {0, 1}, {1, 2}, {2, 3}, {3, 0},
273 // Vertical edges
274 {0, 4}, {1, 5}, {2, 6}, {3, 7},
275 // Top face edges
276 {4, 5}, {5, 6}, {6, 7}, {7, 4},
277 // Face centers
278 {0, 3, 2, 1}, {0, 1, 5, 4},
279 {1, 2, 6, 5}, {2, 3, 7, 6},
280 {0, 4, 7, 3}, {4, 5, 6, 7},
281 // Body center
282 {0, 1, 2, 3, 4, 5, 6, 7}}};
283
284 /// @brief Element type of each elevation span
285 static constexpr std::array<ElemType, 19> elevNodeSpanTypes = {
286 Line2, Line2, Line2, Line2, // Bottom edges
287 Line2, Line2, Line2, Line2, // Vertical edges
288 Line2, Line2, Line2, Line2, // Top edges
289 Quad4, Quad4, Quad4, Quad4, Quad4, Quad4, // Face centers
290 Hex8}; // Body center
291
292 // ============================================================
293 // VTK/Visualization Support
294 // ============================================================
295
296 /// @brief VTK cell type identifier (12 = VTK_HEXAHEDRON)
297 static constexpr int vtkCellType = 12;
298
299 /**
300 * @brief VTK node ordering map
301 *
302 * VTK uses the same ordering as DNDS for Hex8:
303 * VTK nodes 0-7 = corner nodes 0-7
304 */
305 static constexpr std::array<int, 8> vtkNodeOrder = {0, 1, 2, 3, 4, 5, 6, 7};
306 };
307
308} // namespace DNDS::Geom::Elem
Core type aliases, constants, and metaprogramming utilities for the DNDS framework.
#define DNDS_DEVICE_CALLABLE
Definition Defines.hpp:76
int32_t t_index
Definition Geometric.hpp:6
double t_real
Definition Geometric.hpp:8
static constexpr ElemType GetFaceType(t_index)
Get the element type of a face.
Definition Hex8.hpp:229
static DNDS_DEVICE_CALLABLE void Diff1(const TPoint &p, TArray &&v)
Definition Hex8.hpp:48
static DNDS_DEVICE_CALLABLE void Diff0(const TPoint &p, TArray &&v)
Definition Hex8.hpp:22
static DNDS_DEVICE_CALLABLE void Diff3(const TPoint &p, TArray &&v)
Definition Hex8.hpp:147
static DNDS_DEVICE_CALLABLE void Diff2(const TPoint &p, TArray &&v)
Definition Hex8.hpp:100
Eigen::Matrix< real, 5, 1 > v
double order
Definition test_ODE.cpp:257