DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
Line3.hpp
Go to the documentation of this file.
1#pragma once
2// Auto-generated by tools/gen_shape_functions -- DO NOT EDIT
3// Element: Line3
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> struct ShapeFuncImpl;
16
17 // <GEN_SHAPE_FUNCS_BEGIN>
18 template <>
20 {
21 template <class TPoint, class TArray>
22 DNDS_DEVICE_CALLABLE static inline void Diff0(const TPoint &p, TArray &&v)
23 {
24 t_real xi = p[0];
25 const t_real _t0 = ((0.5))*xi;
26 v(0, 0) = _t0*(xi - 1);
27 v(0, 1) = _t0*(xi + 1);
28 v(0, 2) = 1 - ((xi) * (xi));
29 }
30
31 template <class TPoint, class TArray>
32 DNDS_DEVICE_CALLABLE static inline void Diff1(const TPoint &p, TArray &&v)
33 {
34 t_real xi = p[0];
35 v(0, 0) = xi + (-0.5);
36 v(0, 1) = xi + (0.5);
37 v(0, 2) = -2*xi;
38 }
39
40 template <class TPoint, class TArray>
41 DNDS_DEVICE_CALLABLE static inline void Diff2(const TPoint &p, TArray &&v)
42 {
43 t_real xi = p[0];
44 v(0, 0) = 1;
45 v(0, 1) = 1;
46 v(0, 2) = -2;
47 }
48
49 template <class TPoint, class TArray>
50 DNDS_DEVICE_CALLABLE static inline void Diff3(const TPoint &p, TArray &&v)
51 {
52 t_real xi = p[0];
53 // all zero
54 }
55 };
56 // <GEN_SHAPE_FUNCS_END>
57
58
59 /**
60 * @brief Element traits for 3-node quadratic line element (Line3)
61 *
62 * Line3 is a 1D quadratic element with nodes at the endpoints and midpoint.
63 * Used as edge/face elements for 2D/3D high-order meshes.
64 *
65 * Geometry:
66 * - Reference line: xi in [-1, 1]
67 * - Node 0 at xi = -1, Node 1 at xi = +1, Node 2 at xi = 0
68 */
69 template <>
71 {
72 // ============================================================
73 // Core Element Identification
74 // ============================================================
75
76 /// @brief Element type enum value for static dispatch
77 static constexpr ElemType elemType = Line3;
78
79 /// @brief Spatial dimension (1 = line element)
80 static constexpr int dim = 1;
81
82 /// @brief Polynomial order of shape functions (2 = quadratic)
83 static constexpr int order = 2;
84
85 /// @brief Number of vertices (2 endpoints, same as linear)
86 static constexpr int numVertices = 2;
87
88 /// @brief Total number of nodes (3 nodes: 2 vertices + 1 midpoint)
89 static constexpr int numNodes = 3;
90
91 /// @brief Number of faces (0 for 1D element)
92 static constexpr int numFaces = 0;
93
94 /// @brief Parametric space type (line space: xi in [-1, 1])
95 static constexpr ParamSpace paramSpace = LineSpace;
96
97 /// @brief Volume of parametric space (length = 2.0)
98 static constexpr t_real paramSpaceVol = 2.0;
99
100 // ============================================================
101 // Geometry Definition
102 // ============================================================
103
104 /**
105 * @brief Standard coordinates of nodes in parametric space
106 *
107 * Format: {x, y, z} for each node. 3D coordinates for consistency.
108 * Node 0: xi = -1 -> (-1, 0, 0) - left endpoint
109 * Node 1: xi = +1 -> ( 1, 0, 0) - right endpoint
110 * Node 2: xi = 0 -> ( 0, 0, 0) - midpoint
111 */
112 static constexpr std::array<t_real, 3 * 3> standardCoords = {
113 -1, 0, 0, // Node 0: left endpoint
114 1, 0, 0, // Node 1: right endpoint
115 0, 0, 0}; // Node 2: midpoint
116
117 // ============================================================
118 // Face/Edge Type Queries
119 // ============================================================
120
121 /**
122 * @brief Get the element type of a face
123 * @param iFace Face index (unused for Line3 as it has no faces)
124 * @return UnknownElem (1D elements don't have faces)
125 */
126 static constexpr ElemType GetFaceType(t_index /*iFace*/) { return UnknownElem; }
127
128 // ============================================================
129 // Order Elevation (P-Refinement)
130 // ============================================================
131
132 /// @brief Element type after order elevation (O2 has no higher elevation defined)
133 static constexpr ElemType elevatedType = UnknownElem;
134
135 /// @brief Number of additional nodes created during elevation (none for O2)
136 static constexpr int numElevNodes = 0;
137
138 // ============================================================
139 // Bisection (Adaptive Refinement)
140 // ============================================================
141
142 /// @brief Number of sub-elements created when bisecting (2 Line2 elements)
143 static constexpr int numBisect = 2;
144
145 /// @brief Number of bisection variants (only 1 way to bisect a line)
146 static constexpr int numBisectVariants = 1;
147
148 /**
149 * @brief Get the element type of a sub-element after bisection
150 * @param i Sub-element index (0 or 1)
151 * @return Line2 (both sub-elements are linear lines)
152 */
153 static constexpr ElemType GetBisectElemType(t_index /*i*/) { return Line2; }
154
155 /**
156 * @brief Node indices for each sub-element created by bisection
157 *
158 * Bisecting a quadratic line at the midpoint creates 2 linear lines:
159 * Sub-element 0: nodes {0, 2} - left half
160 * Sub-element 1: nodes {2, 1} - right half
161 */
162 static constexpr std::array<tBisectSub, 2> bisectElements = {{
163 {0, 2}, // Left sub-element: nodes 0 to 2
164 {2, 1}}}; // Right sub-element: nodes 2 to 1
165
166 // ============================================================
167 // VTK/Visualization Support
168 // ============================================================
169
170 /// @brief VTK cell type identifier (4 = VTK_QUADRATIC_EDGE)
171 static constexpr int vtkCellType = 4;
172
173 /**
174 * @brief VTK node ordering map
175 *
176 * VTK uses different ordering for quadratic edge:
177 * VTK node 0 = DNDS node 0 (endpoint)
178 * VTK node 1 = DNDS node 2 (midpoint) - note the swap!
179 * VTK node 2 = DNDS node 1 (endpoint)
180 */
181 static constexpr std::array<int, 3> vtkNodeOrder = {0, 2, 1};
182 };
183
184} // 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)
Get the element type of a sub-element after bisection.
Definition Line3.hpp:153
static constexpr ElemType GetFaceType(t_index)
Get the element type of a face.
Definition Line3.hpp:126
static DNDS_DEVICE_CALLABLE void Diff2(const TPoint &p, TArray &&v)
Definition Line3.hpp:41
static DNDS_DEVICE_CALLABLE void Diff1(const TPoint &p, TArray &&v)
Definition Line3.hpp:32
static DNDS_DEVICE_CALLABLE void Diff3(const TPoint &p, TArray &&v)
Definition Line3.hpp:50
static DNDS_DEVICE_CALLABLE void Diff0(const TPoint &p, TArray &&v)
Definition Line3.hpp:22
Eigen::Matrix< real, 5, 1 > v
double order
Definition test_ODE.cpp:257