DNDSR 0.2.1
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
Quad4.hpp
Go to the documentation of this file.
1#pragma once
2// Auto-generated by tools/gen_shape_functions -- DO NOT EDIT
3// Element: Quad4
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 - 1;
28 const t_real _t1 = ((0.25)) * et + (-0.25);
29 const t_real _t2 = xi + 1;
30 const t_real _t3 = ((0.25)) * et + (0.25);
31 v(0, 0) = _t0 * _t1;
32 v(0, 1) = -_t1 * _t2;
33 v(0, 2) = _t2 * _t3;
34 v(0, 3) = -_t0 * _t3;
35 }
36
37 template <class TPoint, class TArray>
38 DNDS_DEVICE_CALLABLE static void Diff1(const TPoint &p, TArray &&v)
39 {
40 t_real xi = p[0];
41 t_real et = p[1];
42 const t_real _t0 = ((0.25)) * et;
43 const t_real _t1 = _t0 + (-0.25);
44 const t_real _t2 = _t0 + (0.25);
45 const t_real _t3 = ((0.25)) * xi;
46 const t_real _t4 = _t3 + (-0.25);
47 const t_real _t5 = _t3 + (0.25);
48 v(0, 0) = _t1;
49 v(0, 1) = -_t1;
50 v(0, 2) = _t2;
51 v(0, 3) = -_t2;
52 v(1, 0) = _t4;
53 v(1, 1) = -_t5;
54 v(1, 2) = _t5;
55 v(1, 3) = -_t4;
56 }
57
58 template <class TPoint, class TArray>
59 DNDS_DEVICE_CALLABLE static void Diff2(const TPoint &p, TArray &&v)
60 {
61 t_real xi = p[0];
62 t_real et = p[1];
63 v(1, 0) = (0.25);
64 v(1, 1) = (-0.25);
65 v(1, 2) = (0.25);
66 v(1, 3) = (-0.25);
67 }
68
69 template <class TPoint, class TArray>
70 DNDS_DEVICE_CALLABLE static void Diff3(const TPoint &p, TArray &&v)
71 {
72 t_real xi = p[0];
73 t_real et = p[1];
74 // all zero
75 }
76 };
77 // <GEN_SHAPE_FUNCS_END>
78
79 // <GEN_ELEM_TRAITS_BEGIN>
80
81 template <>
83 {
84 static constexpr ElemType elemType = Quad4;
85 static constexpr int dim = 2;
86 static constexpr int order = 1;
87 static constexpr int numVertices = 4;
88 static constexpr int numNodes = 4;
89 static constexpr int numFaces = 4;
90 static constexpr ParamSpace paramSpace = QuadSpace;
91 static constexpr t_real paramSpaceVol = 4.0;
92 // 3 * NNodes is a compile-time constant; no overflow possible.
93 // NOLINTNEXTLINE(bugprone-implicit-widening-of-multiplication-result)
94 static constexpr std::array<t_real, 3 * 4> standardCoords = {
95 -1, -1, 0, // Node 0: vertex
96 1, -1, 0, // Node 1: vertex
97 1, 1, 0, // Node 2: vertex
98 -1, 1, 0}; // Node 3: vertex
99
100 static constexpr ElemType GetFaceType(t_index /*iFace*/) { return Line2; }
101
102 static constexpr std::array<std::array<t_index, 10>, 4> faceNodes = {{{0, 1},
103 {1, 2},
104 {2, 3},
105 {3, 0}}};
106
107 static constexpr ElemType elevatedType = Quad9;
108 static constexpr int numElevNodes = 5;
109
110 static constexpr std::array<tElevSpan, 5> elevSpans = {{{0, 1},
111 {1, 2},
112 {2, 3},
113 {3, 0},
114 {0, 1, 2, 3}}};
115
116 static constexpr std::array<ElemType, 5> elevNodeSpanTypes = {
118
119 static constexpr int vtkCellType = 9;
120
121 static constexpr std::array<int, 4> vtkNodeOrder = {0, 1, 2, 3};
122 };
123 // <GEN_ELEM_TRAITS_END>
124
125} // 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 GetFaceType(t_index)
Definition Quad4.hpp:100
static DNDS_DEVICE_CALLABLE void Diff2(const TPoint &p, TArray &&v)
Definition Quad4.hpp:59
static DNDS_DEVICE_CALLABLE void Diff3(const TPoint &p, TArray &&v)
Definition Quad4.hpp:70
static DNDS_DEVICE_CALLABLE void Diff1(const TPoint &p, TArray &&v)
Definition Quad4.hpp:38
static DNDS_DEVICE_CALLABLE void Diff0(const TPoint &p, TArray &&v)
Definition Quad4.hpp:23
Eigen::Matrix< real, 5, 1 > v
double order
Definition test_ODE.cpp:257
const tPoint const tPoint const tPoint & p