Class documentation of Concepts

Loading...
Searching...
No Matches
neumannTraceElement3d.hh
Go to the documentation of this file.
1
7#ifndef hp3dntelement3d_hh
8#define hp3dntelement3d_hh
9
10#include "toolbox/sequence.hh"
12
14#include "space/space.hh"
15#include "space/element.hh"
16#include "hp2D/quad.hh"
17#include "hp3D/hexFunctions.hh"
18#include "hp3D/element.hh"
19
20namespace hp3D{
21
22 using concepts::Real;
23
24 // forward declaration
25 class Hexahedron;
26
27
28 // **************************************************** NeumannTraceElement **
29
41 template<class F = Real>
43 public:
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// */
55// class NTShapeFunction : public concepts::ShapeFunction1D<Real> {
56// public:
57// /** Constructor
58//
59// @param N Number of Shapefunctions
60// @param nP Number of points in which the shape functions are evaluated
61// */
62// NTShapeFunction(const uint N, const uint nP, Real* values)
63// : concepts::ShapeFunction1D<Real>(N, nP)
64// {
65// values_ = values;
66// }
67//
68//
69// protected:
70// virtual std::ostream& info(std::ostream& os) const {
71// return os << "NeumannTraceElement::NTShapeFunction(n = " << n()
72// << ", nP = " << nP() << ", values = "
73// << concepts::Array<Real>(values_, n()*nP()) << ')';
74// }
75// };
76
77
78
79 struct mapPart{
80
81 mapPart(): coarseCell(0){};
82
83 //control constructor: if pointer !=0 work as copy constructor, else like default constructor
84 mapPart(const mapPart* mp): coarseCell((mp)?mp->coarseCell:0){
85 if(mp){
86 t0 = mp->t0;
87 t1 = mp->t1;
88 }
89
90 };
91
92 mapPart(const concepts::Cell2* coarseCellpart, Real t0, Real t1 ):coarseCell(coarseCellpart),t0(t0), t1(t1){};
93
94
95
96
97 //coarse cell mapping
98 const concepts::Cell2* coarseCell;
99 Real t0;
100 Real t1;
101
102 };
103
104
105
106
120 NeumannTraceElement3d( concepts::QuadNd& cell,const ushort* p /*,const mapPart* coarseCell=0*/)
121 : hp2D::Quad<F>(cell,p,0,0), T_(0), shpfct_(nullptr), isComputed_(false){};
122
123 //Decontructor
124 virtual ~NeumannTraceElement3d();
125
126 virtual const concepts::TMatrix<F>& T() const { return T_; }
127
128 /* Uses instead method of parent class
129 virtual const concepts::ElementGraphics<F>* graphics() const{};
130 */
131
140
141
142
143
144
148 void buildT();
149
150// /// Returns the shape functions
151// const concepts::ShapeFunction1D<Real>* shpfct() {
152// //if not yet computed, first time compute them
153// if(!isComputed_)
154// recomputeShapefunctions();
155// return shpfct_.get(); }
156
160 const concepts::Sequence<UnderlyingElement>& uelm()const { return uelm_;}
161
162
170 void addElement(const hp3D::Hexahedron& hex, uint k,
171 Real weight = 1.0);
172
173
174
175 protected:
176 virtual std::ostream& info(std::ostream& os) const;
177
178
179 private:
180
181 // T matrix of the element
183
184 // The shape functions, computed on quadrature points
185 std::unique_ptr<concepts::ShapeFunction1D<Real> > shpfct_;
186
187 //holding information if element's shapefunctions were already evaluated, this is false by default
188 bool isComputed_;
189
191
199 mutable concepts::Array<Real> values_;
200
203
204
205 mapPart coarseCell_;
206 //flag if element has irregular uelms
207 // bool irregular_;
208
216 void compute_(const concepts::Real2d &q, concepts::Array<Real>& val) const;
217
218
219
227 // void irregularCompute_(const concepts::Array<Real>& q, concepts::Array<Real>& val) const;
228
229 };
230
231} // namespace hp2D
232
233#endif // hp2dntelement_hh
234
235
236
Two dimensional cell.
Definition cell.hh:89
virtual const concepts::QuadNd & cell() const
Returns the cell on which the element is built.
Definition quad.hh:196
const ushort * p() const
Definition quad.hh:229
Quad(concepts::QuadNd &cell, const ushort *p, concepts::TColumn< F > *T0, concepts::TColumn< F > *T1)
void computeShpfctVals(const concepts::Real2d &xP, concepts::Array< Real > &val) const
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
NeumannTraceElement3d(concepts::QuadNd &cell, const ushort *p)
virtual const concepts::TMatrix< F > & T() const
Returns the T matrix of the element.
const concepts::Sequence< UnderlyingElement > & uelm() const
void addElement(const hp3D::Hexahedron &hex, uint k, Real weight=1.0)
double Real
Definition typedefs.hh:39
Definition meshDX.hh:23