DNDSR
0.2.1
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
Tet10.hpp
Go to the documentation of this file.
1
#pragma once
2
// Auto-generated by tools/gen_shape_functions -- DO NOT EDIT
3
// Element: Tet10
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
"
9
#include "
Geom/ElementTraitsBase.hpp
"
10
11
namespace
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
<>
20
struct
ShapeFuncImpl
<
Tet10
>
21
{
22
template
<
class
TPo
int
,
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 = et + xi + zt - 1;
29
const
t_real
_t1 = 2 * et;
30
const
t_real
_t2 = 2 * zt;
31
const
t_real
_t3 = 2 * xi - 1;
32
const
t_real
_t4 = -_t0;
33
const
t_real
_t5 = 4 * xi;
34
const
t_real
_t6 = 4 * _t4;
35
v
(0, 0) = _t0 * (_t1 + _t2 + _t3);
36
v
(0, 1) = _t3 * xi;
37
v
(0, 2) = et * (_t1 - 1);
38
v
(0, 3) = zt * (_t2 - 1);
39
v
(0, 4) = _t4 * _t5;
40
v
(0, 5) = _t5 * et;
41
v
(0, 6) = _t6 * et;
42
v
(0, 7) = _t6 * zt;
43
v
(0, 8) = _t5 * zt;
44
v
(0, 9) = 4 * et * zt;
45
}
46
47
template
<
class
TPo
int
,
class
TArray>
48
DNDS_DEVICE_CALLABLE
static
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 = 4 * xi;
54
const
t_real
_t1 = 4 * et;
55
const
t_real
_t2 = 4 * zt;
56
const
t_real
_t3 = _t1 + _t2;
57
const
t_real
_t4 = _t0 + _t3 - 3;
58
const
t_real
_t5 = -_t1;
59
const
t_real
_t6 = -_t2;
60
const
t_real
_t7 = -_t0;
61
const
t_real
_t8 = _t0 - 4;
62
v
(0, 0) = _t4;
63
v
(0, 1) = _t0 - 1;
64
v
(0, 4) = -_t3 - 8 * xi + 4;
65
v
(0, 5) = _t1;
66
v
(0, 6) = _t5;
67
v
(0, 7) = _t6;
68
v
(0, 8) = _t2;
69
v
(1, 0) = _t4;
70
v
(1, 2) = _t1 - 1;
71
v
(1, 4) = _t7;
72
v
(1, 5) = _t0;
73
v
(1, 6) = -_t2 - _t8 - 8 * et;
74
v
(1, 7) = _t6;
75
v
(1, 9) = _t2;
76
v
(2, 0) = _t4;
77
v
(2, 3) = _t2 - 1;
78
v
(2, 4) = _t7;
79
v
(2, 6) = _t5;
80
v
(2, 7) = -_t1 - _t8 - 8 * zt;
81
v
(2, 8) = _t0;
82
v
(2, 9) = _t1;
83
}
84
85
template
<
class
TPo
int
,
class
TArray>
86
DNDS_DEVICE_CALLABLE
static
void
Diff2
(
const
TPoint &
p
, TArray &&
v
)
87
{
88
t_real
xi =
p
[0];
89
t_real
et =
p
[1];
90
t_real
zt =
p
[2];
91
v
(0, 0) = 4;
92
v
(0, 1) = 4;
93
v
(0, 4) = -8;
94
v
(1, 0) = 4;
95
v
(1, 2) = 4;
96
v
(1, 6) = -8;
97
v
(2, 0) = 4;
98
v
(2, 3) = 4;
99
v
(2, 7) = -8;
100
v
(3, 0) = 4;
101
v
(3, 4) = -4;
102
v
(3, 5) = 4;
103
v
(3, 6) = -4;
104
v
(4, 0) = 4;
105
v
(4, 6) = -4;
106
v
(4, 7) = -4;
107
v
(4, 9) = 4;
108
v
(5, 0) = 4;
109
v
(5, 4) = -4;
110
v
(5, 7) = -4;
111
v
(5, 8) = 4;
112
}
113
114
template
<
class
TPo
int
,
class
TArray>
115
DNDS_DEVICE_CALLABLE
static
void
Diff3
(
const
TPoint &
p
, TArray &&
v
)
116
{
117
t_real
xi =
p
[0];
118
t_real
et =
p
[1];
119
t_real
zt =
p
[2];
120
// all zero
121
}
122
};
123
// <GEN_SHAPE_FUNCS_END>
124
125
// <GEN_ELEM_TRAITS_BEGIN>
126
127
template
<>
128
struct
ElementTraits
<
Tet10
>
129
{
130
static
constexpr
ElemType
elemType =
Tet10
;
131
static
constexpr
int
dim = 3;
132
static
constexpr
int
order
= 2;
133
static
constexpr
int
numVertices = 4;
134
static
constexpr
int
numNodes = 10;
135
static
constexpr
int
numFaces = 4;
136
static
constexpr
int
numEdges = 6;
137
static
constexpr
ParamSpace
paramSpace =
TetSpace
;
138
static
constexpr
t_real
paramSpaceVol = 1.0 / 6.0;
139
// 3 * NNodes is a compile-time constant; no overflow possible.
140
// NOLINTNEXTLINE(bugprone-implicit-widening-of-multiplication-result)
141
static
constexpr
std::array<t_real, 3 * 10> standardCoords = {
142
0, 0, 0,
// Node 0: vertex
143
1, 0, 0,
// Node 1: vertex
144
0, 1, 0,
// Node 2: vertex
145
0, 0, 1,
// Node 3: vertex
146
0.5, 0, 0,
// Node 4
147
0.5, 0.5, 0,
// Node 5
148
0, 0.5, 0,
// Node 6
149
0, 0, 0.5,
// Node 7
150
0.5, 0, 0.5,
// Node 8
151
0, 0.5, 0.5};
// Node 9
152
153
static
constexpr
ElemType
GetFaceType
(
t_index
/*iFace*/
) {
return
Tri6
; }
154
155
static
constexpr
std::array<std::array<t_index, 10>, 4> faceNodes = {{{0, 2, 1, 6, 5, 4},
156
{0, 1, 3, 4, 8, 7},
157
{1, 2, 3, 5, 9, 8},
158
{2, 0, 3, 6, 7, 9}}};
159
160
static
constexpr
ElemType
GetEdgeType
(
t_index
/*iEdge*/
) {
return
Line3
; }
161
162
static
constexpr
std::array<std::array<t_index, 3>, 6> edgeNodes = {{{0, 1, 4},
163
{1, 2, 5},
164
{2, 0, 6},
165
{0, 3, 7},
166
{1, 3, 8},
167
{2, 3, 9}}};
168
169
static
constexpr
ElemType
elevatedType =
UnknownElem
;
170
static
constexpr
int
numElevNodes = 0;
171
172
static
constexpr
int
numBisect = 8;
173
static
constexpr
int
numBisectVariants = 3;
174
175
static
constexpr
ElemType
GetBisectElemType
(
t_index
/*i*/
) {
return
Tet4
; }
176
177
static
constexpr
std::array<tBisectSub, 24> bisectElements = {{
// Variant 0
178
{0, 4, 6, 7},
179
{4, 1, 5, 8},
180
{6, 5, 2, 9},
181
{9, 7, 8, 3},
182
{4, 9, 6, 7},
183
{4, 8, 9, 7},
184
{4, 9, 8, 5},
185
{4, 6, 9, 5},
186
// Variant 1
187
{0, 4, 6, 7},
188
{4, 1, 5, 8},
189
{6, 5, 2, 9},
190
{9, 7, 8, 3},
191
{5, 6, 7, 9},
192
{5, 7, 8, 9},
193
{5, 8, 7, 4},
194
{5, 7, 6, 4},
195
// Variant 2
196
{0, 4, 6, 7},
197
{4, 1, 5, 8},
198
{6, 5, 2, 9},
199
{9, 7, 8, 3},
200
{6, 7, 8, 9},
201
{6, 8, 5, 9},
202
{6, 8, 7, 4},
203
{6, 5, 8, 4}}};
204
205
static
constexpr
int
vtkCellType = 24;
206
207
static
constexpr
std::array<int, 10> vtkNodeOrder = {0, 1, 2, 3, 4, 5, 6, 7, 8, 9};
208
};
209
// <GEN_ELEM_TRAITS_END>
210
211
}
// namespace DNDS::Geom::Elem
Defines.hpp
Core type aliases, constants, and metaprogramming utilities for the DNDS framework.
DNDS_DEVICE_CALLABLE
#define DNDS_DEVICE_CALLABLE
Definition
Defines.hpp:76
ElemEnum.hpp
ElementTraitsBase.hpp
Geometric.hpp
DNDS::Geom::Elem
Definition
Elements.cpp:4
DNDS::Geom::Elem::ElemType
ElemType
Definition
ElemEnum.hpp:19
DNDS::Geom::Elem::Tri6
@ Tri6
Definition
ElemEnum.hpp:25
DNDS::Geom::Elem::Tet10
@ Tet10
Definition
ElemEnum.hpp:30
DNDS::Geom::Elem::Line3
@ Line3
Definition
ElemEnum.hpp:22
DNDS::Geom::Elem::UnknownElem
@ UnknownElem
Definition
ElemEnum.hpp:20
DNDS::Geom::Elem::Tet4
@ Tet4
Definition
ElemEnum.hpp:29
DNDS::Geom::Elem::ParamSpace
ParamSpace
Definition
ElemEnum.hpp:42
DNDS::Geom::Elem::TetSpace
@ TetSpace
Definition
ElemEnum.hpp:49
DNDS::Geom::t_index
int32_t t_index
Definition
Geometric.hpp:6
DNDS::Geom::t_real
double t_real
Definition
Geometric.hpp:8
DNDS::Geom::Elem::ElementTraits< Tet10 >::GetEdgeType
static constexpr ElemType GetEdgeType(t_index)
Definition
Tet10.hpp:160
DNDS::Geom::Elem::ElementTraits< Tet10 >::GetBisectElemType
static constexpr ElemType GetBisectElemType(t_index)
Definition
Tet10.hpp:175
DNDS::Geom::Elem::ElementTraits< Tet10 >::GetFaceType
static constexpr ElemType GetFaceType(t_index)
Definition
Tet10.hpp:153
DNDS::Geom::Elem::ElementTraits
Definition
ElementTraitsBase.hpp:22
DNDS::Geom::Elem::ShapeFuncImpl< Tet10 >::Diff2
static DNDS_DEVICE_CALLABLE void Diff2(const TPoint &p, TArray &&v)
Definition
Tet10.hpp:86
DNDS::Geom::Elem::ShapeFuncImpl< Tet10 >::Diff3
static DNDS_DEVICE_CALLABLE void Diff3(const TPoint &p, TArray &&v)
Definition
Tet10.hpp:115
DNDS::Geom::Elem::ShapeFuncImpl< Tet10 >::Diff1
static DNDS_DEVICE_CALLABLE void Diff1(const TPoint &p, TArray &&v)
Definition
Tet10.hpp:48
DNDS::Geom::Elem::ShapeFuncImpl< Tet10 >::Diff0
static DNDS_DEVICE_CALLABLE void Diff0(const TPoint &p, TArray &&v)
Definition
Tet10.hpp:23
DNDS::Geom::Elem::ShapeFuncImpl
Definition
Elements.hpp:174
v
Eigen::Matrix< real, 5, 1 > v
Definition
test_ArrayDOF.cpp:468
order
double order
Definition
test_ODE.cpp:257
p
const tPoint const tPoint const tPoint & p
Definition
test_Reconstruction.cpp:484
src
Geom
Elements
Tet10.hpp
Generated by
1.9.8