DNDSR 0.2.1
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>
16 struct ShapeFuncImpl;
17
18 // <GEN_SHAPE_FUNCS_BEGIN>
19 template <>
21 {
22 template <class TPoint, class TArray>
23 DNDS_DEVICE_CALLABLE static void Diff0(const TPoint &p, TArray &&v)
24 {
25 t_real xi = p[0];
26 t_real et = p[1];
27 t_real zt = p[2];
28 const t_real _t0 = xi - 1;
29 const t_real _t1 = et - 1;
30 const t_real _t2 = ((0.125)) * zt + (-0.125);
31 const t_real _t3 = _t1 * _t2;
32 const t_real _t4 = xi + 1;
33 const t_real _t5 = et + 1;
34 const t_real _t6 = _t2 * _t5;
35 const t_real _t7 = ((0.125)) * zt + (0.125);
36 const t_real _t8 = _t1 * _t7;
37 const t_real _t9 = _t5 * _t7;
38 v(0, 0) = -_t0 * _t3;
39 v(0, 1) = _t3 * _t4;
40 v(0, 2) = -_t4 * _t6;
41 v(0, 3) = _t0 * _t6;
42 v(0, 4) = _t0 * _t8;
43 v(0, 5) = -_t4 * _t8;
44 v(0, 6) = _t4 * _t9;
45 v(0, 7) = -_t0 * _t9;
46 }
47
48 template <class TPoint, class TArray>
49 DNDS_DEVICE_CALLABLE static void Diff1(const TPoint &p, TArray &&v)
50 {
51 t_real xi = p[0];
52 t_real et = p[1];
53 t_real zt = p[2];
54 const t_real _t0 = et - 1;
55 const t_real _t1 = ((0.125)) * zt + (-0.125);
56 const t_real _t2 = _t0 * _t1;
57 const t_real _t3 = et + 1;
58 const t_real _t4 = _t1 * _t3;
59 const t_real _t5 = ((0.125)) * zt + (0.125);
60 const t_real _t6 = _t0 * _t5;
61 const t_real _t7 = _t3 * _t5;
62 const t_real _t8 = xi - 1;
63 const t_real _t9 = _t1 * _t8;
64 const t_real _t10 = xi + 1;
65 const t_real _t11 = _t1 * _t10;
66 const t_real _t12 = _t5 * _t8;
67 const t_real _t13 = _t10 * _t5;
68 const t_real _t14 = ((0.125)) * _t0;
69 const t_real _t15 = _t14 * _t8;
70 const t_real _t16 = _t10 * _t14;
71 const t_real _t17 = ((0.125)) * _t3;
72 const t_real _t18 = _t10 * _t17;
73 const t_real _t19 = _t17 * _t8;
74 v(0, 0) = -_t2;
75 v(0, 1) = _t2;
76 v(0, 2) = -_t4;
77 v(0, 3) = _t4;
78 v(0, 4) = _t6;
79 v(0, 5) = -_t6;
80 v(0, 6) = _t7;
81 v(0, 7) = -_t7;
82 v(1, 0) = -_t9;
83 v(1, 1) = _t11;
84 v(1, 2) = -_t11;
85 v(1, 3) = _t9;
86 v(1, 4) = _t12;
87 v(1, 5) = -_t13;
88 v(1, 6) = _t13;
89 v(1, 7) = -_t12;
90 v(2, 0) = -_t15;
91 v(2, 1) = _t16;
92 v(2, 2) = -_t18;
93 v(2, 3) = _t19;
94 v(2, 4) = _t15;
95 v(2, 5) = -_t16;
96 v(2, 6) = _t18;
97 v(2, 7) = -_t19;
98 }
99
100 template <class TPoint, class TArray>
101 DNDS_DEVICE_CALLABLE static void Diff2(const TPoint &p, TArray &&v)
102 {
103 t_real xi = p[0];
104 t_real et = p[1];
105 t_real zt = p[2];
106 const t_real _t0 = ((0.125)) * zt;
107 const t_real _t1 = _t0 + (-0.125);
108 const t_real _t2 = -_t1;
109 const t_real _t3 = _t0 + (0.125);
110 const t_real _t4 = -_t3;
111 const t_real _t5 = ((0.125)) * xi;
112 const t_real _t6 = _t5 + (-0.125);
113 const t_real _t7 = -_t6;
114 const t_real _t8 = _t5 + (0.125);
115 const t_real _t9 = -_t8;
116 const t_real _t10 = ((0.125)) * et;
117 const t_real _t11 = _t10 + (-0.125);
118 const t_real _t12 = -_t11;
119 const t_real _t13 = _t10 + (0.125);
120 const t_real _t14 = -_t13;
121 v(3, 0) = _t2;
122 v(3, 1) = _t1;
123 v(3, 2) = _t2;
124 v(3, 3) = _t1;
125 v(3, 4) = _t3;
126 v(3, 5) = _t4;
127 v(3, 6) = _t3;
128 v(3, 7) = _t4;
129 v(4, 0) = _t7;
130 v(4, 1) = _t8;
131 v(4, 2) = _t9;
132 v(4, 3) = _t6;
133 v(4, 4) = _t6;
134 v(4, 5) = _t9;
135 v(4, 6) = _t8;
136 v(4, 7) = _t7;
137 v(5, 0) = _t12;
138 v(5, 1) = _t11;
139 v(5, 2) = _t14;
140 v(5, 3) = _t13;
141 v(5, 4) = _t11;
142 v(5, 5) = _t12;
143 v(5, 6) = _t13;
144 v(5, 7) = _t14;
145 }
146
147 template <class TPoint, class TArray>
148 DNDS_DEVICE_CALLABLE static void Diff3(const TPoint &p, TArray &&v)
149 {
150 t_real xi = p[0];
151 t_real et = p[1];
152 t_real zt = p[2];
153 v(9, 0) = (-0.125);
154 v(9, 1) = (0.125);
155 v(9, 2) = (-0.125);
156 v(9, 3) = (0.125);
157 v(9, 4) = (0.125);
158 v(9, 5) = (-0.125);
159 v(9, 6) = (0.125);
160 v(9, 7) = (-0.125);
161 }
162 };
163 // <GEN_SHAPE_FUNCS_END>
164
165 // <GEN_ELEM_TRAITS_BEGIN>
166
167 template <>
169 {
170 static constexpr ElemType elemType = Hex8;
171 static constexpr int dim = 3;
172 static constexpr int order = 1;
173 static constexpr int numVertices = 8;
174 static constexpr int numNodes = 8;
175 static constexpr int numFaces = 6;
176 static constexpr int numEdges = 12;
177 static constexpr ParamSpace paramSpace = HexSpace;
178 static constexpr t_real paramSpaceVol = 8.0;
179 // 3 * NNodes is a compile-time constant; no overflow possible.
180 // NOLINTNEXTLINE(bugprone-implicit-widening-of-multiplication-result)
181 static constexpr std::array<t_real, 3 * 8> standardCoords = {
182 -1, -1, -1, // Node 0: vertex
183 1, -1, -1, // Node 1: vertex
184 1, 1, -1, // Node 2: vertex
185 -1, 1, -1, // Node 3: vertex
186 -1, -1, 1, // Node 4: vertex
187 1, -1, 1, // Node 5: vertex
188 1, 1, 1, // Node 6: vertex
189 -1, 1, 1}; // Node 7: vertex
190
191 static constexpr ElemType GetFaceType(t_index /*iFace*/) { return Quad4; }
192
193 static constexpr std::array<std::array<t_index, 10>, 6> faceNodes = {{{0, 3, 2, 1},
194 {0, 1, 5, 4},
195 {1, 2, 6, 5},
196 {2, 3, 7, 6},
197 {0, 4, 7, 3},
198 {4, 5, 6, 7}}};
199
200 static constexpr ElemType GetEdgeType(t_index /*iEdge*/) { return Line2; }
201
202 static constexpr std::array<std::array<t_index, 3>, 12> edgeNodes = {{{0, 1},
203 {1, 2},
204 {2, 3},
205 {3, 0},
206 {0, 4},
207 {1, 5},
208 {2, 6},
209 {3, 7},
210 {4, 5},
211 {5, 6},
212 {6, 7},
213 {7, 4}}};
214
215 static constexpr ElemType elevatedType = Hex27;
216 static constexpr int numElevNodes = 19;
217
218 static constexpr std::array<tElevSpan, 19> elevSpans = {{{0, 1},
219 {1, 2},
220 {2, 3},
221 {3, 0},
222 {0, 4},
223 {1, 5},
224 {2, 6},
225 {3, 7},
226 {4, 5},
227 {5, 6},
228 {6, 7},
229 {7, 4},
230 {0, 3, 2, 1},
231 {0, 1, 5, 4},
232 {1, 2, 6, 5},
233 {2, 3, 7, 6},
234 {0, 4, 7, 3},
235 {4, 5, 6, 7},
236 {0, 1, 2, 3, 4, 5, 6, 7}}};
237
238 static constexpr std::array<ElemType, 19> elevNodeSpanTypes = {
240
241 static constexpr int vtkCellType = 12;
242
243 static constexpr std::array<int, 8> vtkNodeOrder = {0, 1, 2, 3, 4, 5, 6, 7};
244 };
245 // <GEN_ELEM_TRAITS_END>
246
247} // 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)
Definition Hex8.hpp:191
static constexpr ElemType GetEdgeType(t_index)
Definition Hex8.hpp:200
static DNDS_DEVICE_CALLABLE void Diff1(const TPoint &p, TArray &&v)
Definition Hex8.hpp:49
static DNDS_DEVICE_CALLABLE void Diff0(const TPoint &p, TArray &&v)
Definition Hex8.hpp:23
static DNDS_DEVICE_CALLABLE void Diff3(const TPoint &p, TArray &&v)
Definition Hex8.hpp:148
static DNDS_DEVICE_CALLABLE void Diff2(const TPoint &p, TArray &&v)
Definition Hex8.hpp:101
Eigen::Matrix< real, 5, 1 > v
double order
Definition test_ODE.cpp:257
const tPoint const tPoint const tPoint & p