DNDSR 0.2.1
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
Prism6.hpp
Go to the documentation of this file.
1#pragma once
2// Auto-generated by tools/gen_shape_functions -- DO NOT EDIT
3// Element: Prism6
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 = zt - 1;
29 const t_real _t1 = ((0.5)) * et + ((0.5)) * xi + (-0.5);
30 const t_real _t2 = -(0.5) * _t0;
31 const t_real _t3 = zt + 1;
32 const t_real _t4 = ((0.5)) * _t3;
33 v(0, 0) = _t0 * _t1;
34 v(0, 1) = _t2 * xi;
35 v(0, 2) = _t2 * et;
36 v(0, 3) = -_t1 * _t3;
37 v(0, 4) = _t4 * xi;
38 v(0, 5) = _t4 * et;
39 }
40
41 template <class TPoint, class TArray>
42 DNDS_DEVICE_CALLABLE static void Diff1(const TPoint &p, TArray &&v)
43 {
44 t_real xi = p[0];
45 t_real et = p[1];
46 t_real zt = p[2];
47 const t_real _t0 = ((0.5)) * zt;
48 const t_real _t1 = _t0 + (-0.5);
49 const t_real _t2 = -_t1;
50 const t_real _t3 = _t0 + (0.5);
51 const t_real _t4 = -_t3;
52 const t_real _t5 = ((0.5)) * et;
53 const t_real _t6 = ((0.5)) * xi;
54 const t_real _t7 = _t5 + _t6 + (-0.5);
55 v(0, 0) = _t1;
56 v(0, 1) = _t2;
57 v(0, 3) = _t4;
58 v(0, 4) = _t3;
59 v(1, 0) = _t1;
60 v(1, 2) = _t2;
61 v(1, 3) = _t4;
62 v(1, 5) = _t3;
63 v(2, 0) = _t7;
64 v(2, 1) = -_t6;
65 v(2, 2) = -_t5;
66 v(2, 3) = -_t7;
67 v(2, 4) = _t6;
68 v(2, 5) = _t5;
69 }
70
71 template <class TPoint, class TArray>
72 DNDS_DEVICE_CALLABLE static void Diff2(const TPoint &p, TArray &&v)
73 {
74 t_real xi = p[0];
75 t_real et = p[1];
76 t_real zt = p[2];
77 v(4, 0) = (0.5);
78 v(4, 2) = (-0.5);
79 v(4, 3) = (-0.5);
80 v(4, 5) = (0.5);
81 v(5, 0) = (0.5);
82 v(5, 1) = (-0.5);
83 v(5, 3) = (-0.5);
84 v(5, 4) = (0.5);
85 }
86
87 template <class TPoint, class TArray>
88 DNDS_DEVICE_CALLABLE static void Diff3(const TPoint &p, TArray &&v)
89 {
90 t_real xi = p[0];
91 t_real et = p[1];
92 t_real zt = p[2];
93 // all zero
94 }
95 };
96 // <GEN_SHAPE_FUNCS_END>
97
98 // <GEN_ELEM_TRAITS_BEGIN>
99
100 template <>
102 {
103 static constexpr ElemType elemType = Prism6;
104 static constexpr int dim = 3;
105 static constexpr int order = 1;
106 static constexpr int numVertices = 6;
107 static constexpr int numNodes = 6;
108 static constexpr int numFaces = 5;
109 static constexpr int numEdges = 9;
110 static constexpr ParamSpace paramSpace = PrismSpace;
111 static constexpr t_real paramSpaceVol = 1.0;
112 // 3 * NNodes is a compile-time constant; no overflow possible.
113 // NOLINTNEXTLINE(bugprone-implicit-widening-of-multiplication-result)
114 static constexpr std::array<t_real, 3 * 6> standardCoords = {
115 0, 0, -1, // Node 0: vertex
116 1, 0, -1, // Node 1: vertex
117 0, 1, -1, // Node 2: vertex
118 0, 0, 1, // Node 3: vertex
119 1, 0, 1, // Node 4: vertex
120 0, 1, 1}; // Node 5: vertex
121
122 static constexpr ElemType GetFaceType(t_index iFace)
123 {
124 if (iFace < 3)
125 return Quad4;
126 return Tri3;
127 }
128
129 static constexpr std::array<std::array<t_index, 10>, 5> faceNodes = {{{0, 1, 4, 3},
130 {1, 2, 5, 4},
131 {2, 0, 3, 5},
132 {0, 2, 1},
133 {3, 4, 5}}};
134
135 static constexpr ElemType GetEdgeType(t_index /*iEdge*/) { return Line2; }
136
137 static constexpr std::array<std::array<t_index, 3>, 9> edgeNodes = {{{0, 1},
138 {1, 2},
139 {2, 0},
140 {0, 3},
141 {1, 4},
142 {2, 5},
143 {3, 4},
144 {4, 5},
145 {5, 3}}};
146
147 static constexpr ElemType elevatedType = Prism18;
148 static constexpr int numElevNodes = 12;
149
150 static constexpr std::array<tElevSpan, 12> elevSpans = {{{0, 1},
151 {1, 2},
152 {2, 0},
153 {3, 4},
154 {4, 5},
155 {5, 3},
156 {0, 3},
157 {1, 4},
158 {2, 5},
159 {0, 1, 4, 3},
160 {1, 2, 5, 4},
161 {2, 0, 3, 5}}};
162
163 static constexpr std::array<ElemType, 12> elevNodeSpanTypes = {
165
166 static constexpr int vtkCellType = 13;
167
168 static constexpr std::array<int, 6> vtkNodeOrder = {0, 1, 2, 3, 4, 5};
169 };
170 // <GEN_ELEM_TRAITS_END>
171
172} // 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 GetEdgeType(t_index)
Definition Prism6.hpp:135
static constexpr ElemType GetFaceType(t_index iFace)
Definition Prism6.hpp:122
static DNDS_DEVICE_CALLABLE void Diff1(const TPoint &p, TArray &&v)
Definition Prism6.hpp:42
static DNDS_DEVICE_CALLABLE void Diff2(const TPoint &p, TArray &&v)
Definition Prism6.hpp:72
static DNDS_DEVICE_CALLABLE void Diff3(const TPoint &p, TArray &&v)
Definition Prism6.hpp:88
static DNDS_DEVICE_CALLABLE void Diff0(const TPoint &p, TArray &&v)
Definition Prism6.hpp:23
Eigen::Matrix< real, 5, 1 > v
double order
Definition test_ODE.cpp:257
const tPoint const tPoint const tPoint & p