Class documentation of Concepts

Loading...
Searching...
No Matches
quad.hh
Go to the documentation of this file.
1
6#ifndef quadEdge_hh
7#define quadEdge_hh
8
9#include "basics/typedefs.hh"
12#include "hp2D/quad.hh"
14
15namespace hp2Dedge
16{
17
18 using concepts::Real;
19
20 // ***************************************************** QuadEdgeFunctions **
21
29 {
30 public:
35 QuadEdgeFunctions(const ushort p,
36 const concepts::QuadratureRule2d *intRule);
42 QuadEdgeFunctions(const ushort* p,
43 const concepts::QuadratureRule2d *intRule);
46
49 inline const ushort* p() const
50 {
51 return p_;
52 }
53
54 /* Returns the tangential shape functions in x direction
55 that are multiples of Legendre Polynoms of order 0..p[0]
56 */
57 inline const hp2D::KarniadakisDeriv2* shpfctX_t() const
58 {
59 return shpfctX_t_.get();
60 }
61 /* Returns the tangential shape functions in y direction
62 that are multiples of Legendre Polynoms of order 0..p[1]
63 */
64 inline const hp2D::KarniadakisDeriv2* shpfctY_t() const
65 {
66 return shpfctY_t_.get();
67 }
68 /* Returns the normal shape functions in x direction
69 that are multiples of Integrated Legendre Polynoms of order 1..p[0]
70 order 0 is not included, but not used
71 */
72 inline const concepts::Karniadakis<1, 0>* shpfctX_n() const
73 {
74 return shpfctX_n_.get();
75 }
76 /* Returns the shape functions in y direction
77 that are multiples of Integrated Legendre Polynoms of order 1..p[1]
78 order 0 is not included, but not used
79 */
80 inline const concepts::Karniadakis<1, 0>* shpfctY_n() const
81 {
82 return shpfctY_n_.get();
83 }
84
85 /* Returns the derivatives of the normal shape functions in x direction
86 that are multiples of Legendre Polynoms of order 0..p[0]
87 */
88 inline const concepts::Karniadakis<1, 1>* shpfctDX_n() const
89 {
90 return shpfctDX_n_.get();
91 }
92 /* Returns the derivatives of the normal shape functions in y direction
93 that are multiples of Legendre Polynoms of order 0..p[1]
94 */
95 inline const concepts::Karniadakis<1, 1>* shpfctDY_n() const
96 {
97 return shpfctDY_n_.get();
98 }
99 protected:
102 //void computeShapefunctions_(const concepts::QuadratureRule1d *intX,
103 //const concepts::QuadratureRule1d *intY);
104 private:
106 ushort p_[2];
108 std::unique_ptr<hp2D::KarniadakisDeriv2> shpfctX_t_, shpfctY_t_;
110 std::unique_ptr<concepts::Karniadakis<1, 0> > shpfctX_n_, shpfctY_n_;
112 std::unique_ptr<concepts::Karniadakis<1, 1> > shpfctDX_n_, shpfctDY_n_;
113 };
114
115 // ****************************************************************** Quad **
116
127 template<class F = Real>
128 class Quad : public hp2D::BaseQuad<F>, public QuadEdgeFunctions
129 {
130 public:
139
145 void recomputeShapefunctions(const uint nq[2]);
146 virtual ~Quad();
147 protected:
148 virtual std::ostream& info(std::ostream& os) const;
149 private:
151 static std::unique_ptr<concepts::ElementGraphics<F> > graphics_;
152 };
153
154 // ****************************************************************** Edge **
155
161 template<class F = Real>
162 class Edge : public concepts::Element<F>
163 {
164 public:
169 Edge(const Quad<F>& elm, const ushort k);
170
171 virtual ~Edge();
172
173 inline concepts::Real2d vertex(uint i) const
174 {
176 return i ? chi(1.0) : chi(0.0);
177 }
178
179 virtual const concepts::Edge& support() const
180 {
181 return *elm_.cell().connector().edge(k_);
182 }
183
185 inline const ushort edge() const
186 {
187 return k_;
188 }
189
191 inline const Quad<F>& elm() const
192 {
193 return elm_;
194 }
195
197 inline const ushort direction() const
198 {
199 return l_;
200 }
201
203 inline const Real sign() const
204 {
205 return k_ % 3 == 0 ? -1.0 : 1.0;
206 }
207
208 /* Returns the shape functions
209 that are multiples of Legendre Polynoms of order 0..p
210 */
211 inline const hp2D::KarniadakisDeriv2* shpfct() const
212 {
213 return shpfct_;
214 }
217 {
218 return int_;
219 }
220
222 inline concepts::Real2d localCoords(const Real t) const
223 {
225 x[l_] = x_;
226 x[(l_ + 1) % 2] = t;
227 return x;
228 }
229
232 inline concepts::Real2d chi(const Real t) const
233 {
235 return elm_.chi(x[0], x[1]);
236 }
237
241 inline concepts::MapReal2d jacobian(const Real t) const
242 {
244 return elm_.jacobian(x[0], x[1]);
245 }
246
248 inline concepts::MapReal2d jacobianInverse(const Real t) const
249 {
250 return jacobian(t).inverse();
251 }
252
254 inline Real jacobianDeterminant(const Real t) const
255 {
256 return jacobian(t).determinant();
257 }
258
260 inline Real diffElement(const Real t) const
261 {
262 // factor comes due to integration over [-1,1], but differential element
263 // was defined for reference element [0,1]
264 return elm_.cell().lineElement(t, k_) * 0.5;
265 }
266
268 virtual const concepts::TMatrixBase<F>& T() const
269 {
270 return *T_;
271 }
272 protected:
273 virtual std::ostream& info(std::ostream& os) const;
274 private:
276 const Quad<F>& elm_;
278 const ushort k_;
280 const ushort l_;
286 const hp2D::KarniadakisDeriv2* shpfct_;
288 const concepts::QuadratureRule1d *int_;
290 const concepts::TMatrixBase<F>* T_;
291 };
292
293} // namespace hp2Dedge
294
295#endif // quadEdge_hh
Mapping< F, DimX, DimY > inverse() const
Returns the inverse of the matrix.
F determinant() const
Returns the determinant of the matrix (only valid for square matrices)
virtual const concepts::QuadNd & cell() const
Returns the cell on which the element is built.
Definition quad.hh:196
virtual const concepts::TMatrixBase< F > & T() const
T-Matrix of the appropiate Quad, not used.
Definition quad.hh:268
const ushort edge() const
Returns number of the edge.
Definition quad.hh:185
const Quad< F > & elm() const
Returns element.
Definition quad.hh:191
const Real sign() const
Returns sign of outer normal vector, e.g. left edge -1, right edge +1.
Definition quad.hh:203
concepts::MapReal2d jacobian(const Real t) const
Definition quad.hh:241
const concepts::QuadratureRule1d * integration() const
Returns the integration rule.
Definition quad.hh:216
concepts::Real2d localCoords(const Real t) const
coordinate of point on the edge inside reference element [0,1]^2
Definition quad.hh:222
concepts::Real2d chi(const Real t) const
Definition quad.hh:232
concepts::MapReal2d jacobianInverse(const Real t) const
Computes the inverse of the Jacobian.
Definition quad.hh:248
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
const ushort direction() const
Returns direction of edge on reference quad [0,1]^2, 0 - x, 1 - y.
Definition quad.hh:197
Real diffElement(const Real t) const
Computes the differential element for integration over [-1,1].
Definition quad.hh:260
Edge(const Quad< F > &elm, const ushort k)
Real jacobianDeterminant(const Real t) const
Computes the determinant of the Jacobian.
Definition quad.hh:254
void computeShapefunctions_(const concepts::QuadratureRule2d *intRule)
gets the shapefunctions, used in both constructors
QuadEdgeFunctions(const ushort *p, const concepts::QuadratureRule2d *intRule)
virtual ~QuadEdgeFunctions()
Destructor.
const ushort * p() const
Definition quad.hh:49
QuadEdgeFunctions(const ushort p, const concepts::QuadratureRule2d *intRule)
Quad(concepts::Quad2d &cell, ushort *p, concepts::TColumn< F > *T0, concepts::TColumn< F > *T1)
void recomputeShapefunctions()
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
virtual const concepts::ElementGraphics< F > * graphics() const
Returns element graphics class.
#define conceptsAssert(cond, exc)
double Real
Definition typedefs.hh:39