DNDSR 0.2.1
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
Pyramid5.hpp
Go to the documentation of this file.
1#pragma once
2// Auto-generated by tools/gen_shape_functions -- DO NOT EDIT
3// Element: Pyramid5
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 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 // singularity guard
29 if (std::abs(1 - zt) < 1e-15)
30 zt -= DNDS::signP(1 - zt) * 1e-15;
31 const t_real _t0 = zt - 1;
32 const t_real _t1 = _t0 + et;
33 const t_real _t2 = _t0 + xi;
34 const t_real _t3 = (1.0 / (4 * zt - 4));
35 const t_real _t4 = 1 - zt;
36 const t_real _t5 = _t4 + xi;
37 const t_real _t6 = ((0.25)) / _t0;
38 const t_real _t7 = _t4 + et;
39 v(0, 0) = -_t1 * _t2 * _t3;
40 v(0, 1) = _t1 * _t5 * _t6;
41 v(0, 2) = -_t3 * _t5 * _t7;
42 v(0, 3) = _t2 * _t6 * _t7;
43 v(0, 4) = zt;
44 }
45
46 template <class TPoint, class TArray>
47 static void Diff1(const TPoint &p, TArray &&v)
48 {
49 t_real xi = p[0];
50 t_real et = p[1];
51 t_real zt = p[2];
52 // singularity guard
53 if (std::abs(1 - zt) < 1e-15)
54 zt -= DNDS::signP(1 - zt) * 1e-15;
55 const t_real _t0 = zt - 1;
56 const t_real _t1 = _t0 + et;
57 const t_real _t2 = ((0.25)) / _t0;
58 const t_real _t3 = 1 - zt;
59 const t_real _t4 = _t3 + et;
60 const t_real _t5 = _t0 + xi;
61 const t_real _t6 = _t3 + xi;
62 const t_real _t7 = 2 * zt;
63 const t_real _t8 = et * xi;
64 const t_real _t9 = ((zt) * (zt));
65 const t_real _t10 = -_t7 + _t9 + 1;
66 const t_real _t11 = ((0.25)) / _t10;
67 const t_real _t12 = _t11 * (_t7 + _t8 - _t9 - 1);
68 const t_real _t13 = _t11 * (-_t10 - _t8);
69 v(0, 0) = -_t1 * _t2;
70 v(0, 1) = _t1 * _t2;
71 v(0, 2) = -_t2 * _t4;
72 v(0, 3) = _t2 * _t4;
73 v(1, 0) = -_t2 * _t5;
74 v(1, 1) = _t2 * _t6;
75 v(1, 2) = -_t2 * _t6;
76 v(1, 3) = _t2 * _t5;
77 v(2, 0) = _t12;
78 v(2, 1) = _t13;
79 v(2, 2) = _t12;
80 v(2, 3) = _t13;
81 v(2, 4) = 1;
82 }
83
84 template <class TPoint, class TArray>
85 static void Diff2(const TPoint &p, TArray &&v)
86 {
87 t_real xi = p[0];
88 t_real et = p[1];
89 t_real zt = p[2];
90 // singularity guard
91 if (std::abs(1 - zt) < 1e-15)
92 zt -= DNDS::signP(1 - zt) * 1e-15;
93 const t_real _t0 = ((zt) * (zt));
94 const t_real _t1 = ((zt) * (zt) * (zt));
95 const t_real _t2 = et * xi;
96 const t_real _t3 = -_t2 / (-6 * _t0 + 2 * _t1 + 6 * zt - 2);
97 const t_real _t4 = ((0.5)) * _t2 / (-3 * _t0 + _t1 + 3 * zt - 1);
98 const t_real _t5 = -1 / (4 * zt - 4);
99 const t_real _t6 = zt - 1;
100 const t_real _t7 = ((0.25)) / _t6;
101 const t_real _t8 = ((0.25)) / ((_t6) * (_t6));
102 const t_real _t9 = _t8 * xi;
103 const t_real _t10 = -_t9;
104 const t_real _t11 = _t8 * et;
105 const t_real _t12 = -_t11;
106 v(2, 0) = _t3;
107 v(2, 1) = _t4;
108 v(2, 2) = _t3;
109 v(2, 3) = _t4;
110 v(3, 0) = _t5;
111 v(3, 1) = _t7;
112 v(3, 2) = _t5;
113 v(3, 3) = _t7;
114 v(4, 0) = _t9;
115 v(4, 1) = _t10;
116 v(4, 2) = _t9;
117 v(4, 3) = _t10;
118 v(5, 0) = _t11;
119 v(5, 1) = _t12;
120 v(5, 2) = _t11;
121 v(5, 3) = _t12;
122 }
123
124 template <class TPoint, class TArray>
125 static void Diff3(const TPoint &p, TArray &&v)
126 {
127 t_real xi = p[0];
128 t_real et = p[1];
129 t_real zt = p[2];
130 // singularity guard
131 if (std::abs(1 - zt) < 1e-15)
132 zt -= DNDS::signP(1 - zt) * 1e-15;
133 const t_real _t0 = ((zt) * (zt) * (zt) * (zt));
134 const t_real _t1 = ((zt) * (zt));
135 const t_real _t2 = ((zt) * (zt) * (zt));
136 const t_real _t3 = et * xi;
137 const t_real _t4 = ((1.5)) * _t3 / (_t0 + 6 * _t1 - 4 * _t2 - 4 * zt + 1);
138 const t_real _t5 = -3 * _t3 / (2 * _t0 + 12 * _t1 - 8 * _t2 - 8 * zt + 2);
139 const t_real _t6 = zt - 1;
140 const t_real _t7 = ((0.5)) / ((_t6) * (_t6) * (_t6));
141 const t_real _t8 = _t7 * xi;
142 const t_real _t9 = -_t8;
143 const t_real _t10 = _t7 * et;
144 const t_real _t11 = -_t10;
145 const t_real _t12 = ((0.25)) / ((_t6) * (_t6));
146 const t_real _t13 = -_t12;
147 v(2, 0) = _t4;
148 v(2, 1) = _t5;
149 v(2, 2) = _t4;
150 v(2, 3) = _t5;
151 v(6, 0) = _t9;
152 v(6, 1) = _t8;
153 v(6, 2) = _t9;
154 v(6, 3) = _t8;
155 v(7, 0) = _t11;
156 v(7, 1) = _t10;
157 v(7, 2) = _t11;
158 v(7, 3) = _t10;
159 v(9, 0) = _t12;
160 v(9, 1) = _t13;
161 v(9, 2) = _t12;
162 v(9, 3) = _t13;
163 }
164 };
165 // <GEN_SHAPE_FUNCS_END>
166
167 // <GEN_ELEM_TRAITS_BEGIN>
168
169 template <>
171 {
172 static constexpr ElemType elemType = Pyramid5;
173 static constexpr int dim = 3;
174 static constexpr int order = 1;
175 static constexpr int numVertices = 5;
176 static constexpr int numNodes = 5;
177 static constexpr int numFaces = 5;
178 static constexpr int numEdges = 8;
179 static constexpr ParamSpace paramSpace = PyramidSpace;
180 static constexpr t_real paramSpaceVol = 4.0 / 3.0;
181 // 3 * NNodes is a compile-time constant; no overflow possible.
182 // NOLINTNEXTLINE(bugprone-implicit-widening-of-multiplication-result)
183 static constexpr std::array<t_real, 3 * 5> standardCoords = {
184 -1, -1, 0, // Node 0: vertex
185 1, -1, 0, // Node 1: vertex
186 1, 1, 0, // Node 2: vertex
187 -1, 1, 0, // Node 3: vertex
188 0, 0, 1}; // Node 4: vertex
189
190 static constexpr ElemType GetFaceType(t_index iFace)
191 {
192 if (iFace < 1)
193 return Quad4;
194 return Tri3;
195 }
196
197 static constexpr std::array<std::array<t_index, 10>, 5> faceNodes = {{{0, 3, 2, 1},
198 {0, 1, 4},
199 {1, 2, 4},
200 {2, 3, 4},
201 {3, 0, 4}}};
202
203 static constexpr ElemType GetEdgeType(t_index /*iEdge*/) { return Line2; }
204
205 static constexpr std::array<std::array<t_index, 3>, 8> edgeNodes = {{{0, 1},
206 {1, 2},
207 {2, 3},
208 {3, 0},
209 {0, 4},
210 {1, 4},
211 {2, 4},
212 {3, 4}}};
213
214 static constexpr ElemType elevatedType = Pyramid14;
215 static constexpr int numElevNodes = 9;
216
217 static constexpr std::array<tElevSpan, 9> elevSpans = {{{0, 1},
218 {1, 2},
219 {2, 3},
220 {3, 0},
221 {0, 4},
222 {1, 4},
223 {2, 4},
224 {3, 4},
225 {0, 3, 2, 1}}};
226
227 static constexpr std::array<ElemType, 9> elevNodeSpanTypes = {
229
230 static constexpr int vtkCellType = 14;
231
232 static constexpr std::array<int, 5> vtkNodeOrder = {0, 1, 2, 3, 4};
233 };
234 // <GEN_ELEM_TRAITS_END>
235
236} // namespace DNDS::Geom::Elem
Core type aliases, constants, and metaprogramming utilities for the DNDS framework.
int32_t t_index
Definition Geometric.hpp:6
double t_real
Definition Geometric.hpp:8
constexpr real signP(real a)
"Signum, biased toward +1": treats 0 as positive.
Definition Defines.hpp:549
static constexpr ElemType GetEdgeType(t_index)
Definition Pyramid5.hpp:203
static constexpr ElemType GetFaceType(t_index iFace)
Definition Pyramid5.hpp:190
static void Diff1(const TPoint &p, TArray &&v)
Definition Pyramid5.hpp:47
static void Diff2(const TPoint &p, TArray &&v)
Definition Pyramid5.hpp:85
static void Diff0(const TPoint &p, TArray &&v)
Definition Pyramid5.hpp:23
static void Diff3(const TPoint &p, TArray &&v)
Definition Pyramid5.hpp:125
Eigen::Matrix< real, 5, 1 > v
double order
Definition test_ODE.cpp:257
const tPoint const tPoint const tPoint & p