Class documentation of Concepts

Loading...
Searching...
No Matches
element.hh
Go to the documentation of this file.
1
6#ifndef Element1D_h
7#define Element1D_h
8
9#include "basics/debug.hh"
10#include "basics/output.hh"
11#include "basics/typedefs.hh"
13#include "geometry/connector.hh"
14#include "geometry/cell1D.hh"
15#include "geometry/integral.hh"
16#include "space/element.hh"
17#include "space/tmatrix.hh"
21#include "lineGraphics.hh"
22
23// debugging Element
24#define ElementConstr_D 0
25#define ElementDestr_D 0
26
27namespace hp1D {
28using concepts::Real;
29
30// ********************************************************* IntegrableElm **
31
35public:
37
40 const concepts::Real3d chi(const Real x) const {
41 return cell_.elemMap(x);
42 }
43
45 inline Real jacobianDeterminant(const Real x) const {
46 // factor comes due to integration over [-1,1]^2
47 return cell_.jacobianDeterminant(x) * 0.5;
48 }
49
52 return int_.get();
53 }
54
63 return rule_;
64 }
65
66 virtual bool quadraturePoint(uint i, intPoint& p, intFormType form = ZERO,
67 bool localCoord = false) const;
68protected:
72 std::unique_ptr<concepts::QuadratureRule1d> int_;
73 // Factory for the integration rule
74 static concepts::QuadRuleFactory rule_;
75};
76
77// ********************************************************* BaseElement **
78
79template<class F>
81public:
82 typedef F FieldT;
83
84 BaseElement(const concepts::EdgeNd& cell, ushort p) :
85 IntegrableElm(cell), p_(p) {
86 }
87 ;
88
89 virtual ~BaseElement();
90
91 virtual const concepts::Edge& support() const {
92 return cell_.connector();
93 }
94
95 virtual concepts::Real3d vertex(uint i) const {
96 return cell_.vertex(i);
97 }
98
99 virtual const concepts::EdgeNd& cell() const {
100 return cell_;
101 }
102
103 ushort p() const {
104 return p_;
105 }
106
109
113 virtual void recomputeShapefunctions() = 0;
114
116 virtual const concepts::ShapeFunction1D<Real>* shpfct() const = 0;
117
118protected:
119 virtual std::ostream& info(std::ostream& os) const;
120
122 // concepts::TMatrix<Real> T_;
123
124 ushort p_;
125
126 static std::unique_ptr<LineGraphics> graphics_;
127
128};
129
131
132// *************************************************************** KarniadakisMixin ***
133
136template<class F>
138public:
139 KarniadakisMixin(const concepts::EdgeNd& cell, ushort p) :
140 BaseElement<F>(cell, p), shpfct_(nullptr), shpfctD_(nullptr), shpfctDD_(
141 nullptr), secondDerivativeUsed_(false) {
143 }
144
148 void recomputeShapefunctions() override;
149
151 const concepts::Karniadakis<1, 0>* shpfct() const override {
152 return shpfct_.get();
153 }
156 return shpfctD_.get();
157 }
158
161 if (!secondDerivativeUsed_) {
162 shpfctDD_.reset(
163 new concepts::Karniadakis<1, 2>(this->p_, this->int_->abscissas(),
164 this->int_->n()));
165 secondDerivativeUsed_ = true;
166 }
167 return shpfctDD_.get();
168 }
169
170private:
175 void recomputeSecondDerivativeOfShapefunctions_();
176
178 std::unique_ptr<concepts::Karniadakis<1, 0> > shpfct_;
180 std::unique_ptr<concepts::Karniadakis<1, 1> > shpfctD_;
181
183 mutable std::unique_ptr<concepts::Karniadakis<1, 2> > shpfctDD_;
184
185 mutable bool secondDerivativeUsed_;
186
187};
188
189
190
191// ******************************************************** LegendreMixin ***
192
195template<class F>
196class LegendreMixin: public BaseElement<F> {
197public:
198 LegendreMixin(const concepts::EdgeNd& cell, ushort p) : BaseElement<F>(cell, p) {
200 }
201
205 void recomputeShapefunctions() override;
206
208 const concepts::Legendre* shpfct() const override {
209 return shpfct_.get();
210 }
211
212private:
214 std::unique_ptr<concepts::Legendre> shpfct_;
215
216};
217
218
219// ******************************************************* GenericElement **
220
224template<class BaseT>
225class GenericElement: public BaseT {
226public:
233 GenericElement(const concepts::EdgeNd& cell, uint p,
236 : BaseT(cell,p), T_(T1)
237 {
238 if (T0)
239 T_.append(T0);
240 DEBUGL(ElementConstr_D, cell << ", p = " << p << ", T = " << T_);
241 }
242
243 ~GenericElement() override
244 {
245 DEBUGL(ElementDestr_D, "done.");
246 }
247
248 const concepts::TMatrix<typename BaseT::FieldT>& T() const override {
249 return T_;
250 }
251
256
257protected:
258 virtual std::ostream& info(std::ostream& os) const override;
259
262
263};
264
266template<class F>
268
269
270template<class F>
272
273}// namespace hp1D
274
275namespace concepts {
276
277// ****************************************************************** Scan **
278
280
281template<class F>
282class Scan<hp1D::BaseElement<F> > : public Scan<ElementWithCell<F> > {
283public:
285};
286
287} // namespace concepts
288
289#endif // Element1D_h
virtual Real3d elemMap(const Real coord_local) const
Element map from point local coordinates in 1D.
virtual Real jacobianDeterminant(const Real x) const =0
Returns the determinant of the Jacobian.
Edge & connector() const
Returns the connector (topology)
Definition cell1D.hh:59
virtual Real3d vertex(const uint i) const
Returns the coordinates of the ith vertex.
hp1D::BaseElement< F > & operator++(int)=0
Returns the next element in the scanned set.
void append(TColumn< F > *T)
ushort p_
T matrix of the element.
Definition element.hh:124
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
virtual const concepts::ElementGraphics< Real > * graphics() const
Returns element graphics class.
virtual void recomputeShapefunctions()=0
virtual const concepts::EdgeNd & cell() const
Returns the cell on which the element is built.
Definition element.hh:99
virtual const concepts::ShapeFunction1D< Real > * shpfct() const =0
Returns the shape functions.
concepts::TMatrix< typename BaseT::FieldT > T_
T matrix of the element.
Definition element.hh:261
void appendT(concepts::TColumn< typename BaseT::FieldT > *T)
Appends the T columns to the T matrix.
Definition element.hh:253
GenericElement(const concepts::EdgeNd &cell, uint p, concepts::TColumn< typename BaseT::FieldT > *T0, concepts::TColumn< typename BaseT::FieldT > *T1)
Definition element.hh:233
Real jacobianDeterminant(const Real x) const
Computes the determinant of the Jacobian.
Definition element.hh:45
std::unique_ptr< concepts::QuadratureRule1d > int_
The integration rule.
Definition element.hh:72
const concepts::QuadratureRule1d * integration() const
Returns the integration rule.
Definition element.hh:51
static concepts::QuadRuleFactory & rule()
Definition element.hh:62
const concepts::Real3d chi(const Real x) const
Definition element.hh:40
virtual bool quadraturePoint(uint i, intPoint &p, intFormType form=ZERO, bool localCoord=false) const
const concepts::EdgeNd & cell_
The cell.
Definition element.hh:70
The following two types are shape function mixins.
Definition element.hh:137
void recomputeShapefunctions() override
const concepts::Karniadakis< 1, 1 > * shpfctD() const
Returns the derivatives of the shape functions.
Definition element.hh:155
const concepts::Karniadakis< 1, 2 > * shpfctDD() const
Returns the second derivatives of the shape functions.
Definition element.hh:160
const concepts::Karniadakis< 1, 0 > * shpfct() const override
Returns the shape functions.
Definition element.hh:151
void recomputeShapefunctions() override
const concepts::Legendre * shpfct() const override
Returns the shape functions.
Definition element.hh:208
#define DEBUGL(doit, msg)
Definition debug.hh:40
double Real
Definition typedefs.hh:39