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