Class documentation of Concepts

Loading...
Searching...
No Matches
neumannTraceElement.hh
Go to the documentation of this file.
1
8#ifndef hp2dntelement_hh
9#define hp2dntelement_hh
10
11#include "toolbox/sequence.hh"
13
15#include "space/space.hh"
16#include "space/element.hh"
17#include "hp1D/element.hh"
18#include "hp2D/quad.hh"
19
20namespace hp2D {
21
22using concepts::Real;
23
24// forward declaration
25template<typename F>
26class Quad;
27
28// **************************************************** NeumannTraceElement **
29
41template<class F = Real>
43public:
45
46 /*
47 * Shapefunction as normal derivatives of origin basisfunctions
48 * these class only represents the Shapefunctions evaluated on quadrature Points
49 * since it controls if elements shapefunctions quadrature evaluation is done already
50 * and if not that will be done.
51 *
52 * If this shapefunctions should be evaluated on certain points (e.g. in functions like hp2D::Value)
53 * use the computeShpfct routine
54 */
56 public:
62 NTShapeFunction(const uint N, const uint nP, Real* values) :
63 concepts::ShapeFunction1D<Real>(N, nP) {
65 }
66
67 protected:
68 virtual std::ostream& info(std::ostream& os) const {
69 return os << concepts::typeOf(*this) << "(n = " << n() << ", nP = "
70 << nP() << ", values = " << concepts::Array<Real>(values_, n() * nP())
71 << ')';
72 }
73 };
74
75 //helper structure for coarse to fine edge orientation
76 struct mapPart {
77 //default constructor
78 mapPart() :
79 coarseCell(0) {
80 }
81 ;
82 //control constructor: if pointer !=0 work as copy constructor, else like default constructor
83 mapPart(const mapPart* mp) :
84 coarseCell((mp) ? mp->coarseCell : 0) {
85 if (mp) {
86 t0 = mp->t0;
87 t1 = mp->t1;
88 }
89 }
90
91 mapPart(const concepts::EdgeNd* coarseCellpart, Real t0, Real t1) :
92 coarseCell(coarseCellpart), t0(t0), t1(t1) {
93 }
94 ;
95
96 //coarse cell mapping
97 const concepts::EdgeNd* coarseCell;
98 Real t0;
99 Real t1;
100 };
101
116 const mapPart* coarseCell = 0);
117
118 //Decontructor
119 virtual ~NeumannTraceElement();
120
121 virtual const concepts::TMatrix<F>& T() const {
122 return T_;
123 }
124
133 concepts::Array<Real>& val) const;
134
139
143 void buildT();
144
147
152 return uelm_;
153 }
154
162 void addElement(const hp2D::Quad<Real>& quad, uint k, Real weight = 1.0);
163
170 void addIrrElement(const hp2D::Quad<Real>& coarseQuad, uint k_coarse,
171 Real weight);
172
173protected:
174 virtual std::ostream& info(std::ostream& os) const;
175
176private:
177
178 // The shape functions, computed on quadrature points on demand
179 mutable std::unique_ptr<concepts::ShapeFunction1D<Real> > shpfct_;
180
181 // T matrix of the element
183
184 //holding information if element's shapefunctions were already evaluated, this is false by default
185 //mutable bool isComputed_;
186
187 //flag if element has irregular uelms
188 bool irregular_;
189
190 mapPart coarseCell_;
191
193
201 mutable concepts::Array<Real> values_;
202
205
213 void compute_(const concepts::Array<Real>& q,
214 concepts::Array<Real>& val) const;
215
223 void irregularCompute_(const concepts::Array<Real>& q,
224 concepts::Array<Real>& val) const;
225
226};
227
228} // namespace hp2D
229
230#endif // hp2dntelement_hh
const Real * values() const
Returns the values of the shape functions.
uint n() const
Returns the number of shape functions.
Real * values_
Values of the shape functions.
ShapeFunction1D(const uint n, const uint nP)
virtual const concepts::EdgeNd & cell() const
Returns the cell on which the element is built.
Definition element.hh:99
NTShapeFunction(const uint N, const uint nP, Real *values)
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
virtual const concepts::ShapeFunction1D< Real > * shpfct() const
Returns the shape functions.
void addElement(const hp2D::Quad< Real > &quad, uint k, Real weight=1.0)
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
NeumannTraceElement(const concepts::EdgeNd &cell, uint p, const mapPart *coarseCell=0)
void addIrrElement(const hp2D::Quad< Real > &coarseQuad, uint k_coarse, Real weight)
virtual void recomputeShapefunctions()
virtual const concepts::TMatrix< F > & T() const
Returns the T matrix of the element.
void computeShpfctVals(const concepts::Array< Real > &xP, concepts::Array< Real > &val) const
const concepts::Sequence< UnderlyingElement > & uelm() const
std::string typeOf(const T &t)
Definition output.hh:43
double Real
Definition typedefs.hh:39