Class documentation of Concepts

Loading...
Searching...
No Matches
neumannTrace.hh
Go to the documentation of this file.
1
9#ifndef hp2dneumanntrace_hh
10#define hp2dneumanntrace_hh
11
12#include "toolbox/sequence.hh"
14#include "space/space.hh"
15#include "space/element.hh"
17#include "hp1D/element.hh"
18#include "hp2D/quad.hh"
19
20namespace hp2D {
21
22 using concepts::Real;
23
24 // forward declaration
25 template<typename F>
26 class Quad;
27
28 // **************************************************** NTElement_BA **
29
38 template<class F = Real>
40 , public hp1D::IntegrableElm {
41 public:
43
45 public:
51 ShapeFunction(const uint n, const uint nP, Real* values)
52 : concepts::ShapeFunction1D<Real>(n, nP)
53 {
55 }
56 protected:
57 virtual std::ostream& info(std::ostream& os) const {
58 return os << concepts::typeOf(*this)<<"(n = " << n()
59 << ", nP = " << nP() << ", values = "
60 << concepts::Array<Real>(values_, n()*nP()) << ')';
61 }
62 };
63
69
70
71 virtual ~NTElement_BA();
72
73 virtual const concepts::Edge& support() const {
74 return cell_.connector(); }
75 virtual concepts::Real3d vertex(uint i) const;
76 virtual const concepts::Edge2d& cell() const { return cell_; }
77
78 virtual const concepts::TMatrix<F>& T() const { return T_; }
79
84
87 return shpfct_.get(); }
93 throw(concepts::MissingFeature("not implemented"));
94 return shpfctD_.get();
95 }
96
97
101 const concepts::Sequence<UnderlyingElement>& uelm()const { return uelm_;}
102
113 const concepts::Z2 mapRho() const {return mapRho_;}
114
115
119 void addElement(const hp2D::Quad<Real>& quad, uint k,
120 Real weight = 1.0);
121
124 T_.append(T);
125 }
126
127 ushort p() const {
128 assert(p_.size()>0);
129 if(p_.size()>1)
130 assert(p_[0]==p_[1]);
131 return p_[0];
132 }
133
134 protected:
135 virtual std::ostream& info(std::ostream& os) const;
136
138 //vectorial:: //als unique_ptr durch new routine in tmatrix.hh
140 private:
142 const concepts::Edge2d& cell_;
143
145
147 std::unique_ptr<concepts::ShapeFunction1D<Real> > shpfct_;
149 std::unique_ptr<concepts::ShapeFunction1D<Real> > shpfctD_;
153 concepts::Array<Real> values_, valuesD_;
159 concepts::Sequence<Real> weights_;
160
161 //Orientation of the mapping, the mapping is not oriented like to topological one in general
162 concepts::Z2 mapRho_;
163 };
164
165} // namespace hp2D
166
167#endif // hp2dneumanntrace_hh
Edge & connector() const
Returns the connector (topology)
Definition cell1D.hh:59
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)
void append(TColumn< F > *T)
ShapeFunction(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::TMatrix< F > & T() const
Returns the T matrix of the element.
concepts::TMatrix< Real > T_
T matrix of the element.
const concepts::ShapeFunction1D< Real > * shpfctD() const
virtual const concepts::Edge2d & cell() const
Returns the cell on which the element is built.
void appendT(concepts::TColumn< F > *T)
Appends the T columns to the T matrix.
const concepts::ShapeFunction1D< Real > * shpfct() const
Returns the shape functions.
void addElement(const hp2D::Quad< Real > &quad, uint k, Real weight=1.0)
NTElement_BA(const concepts::Edge2d &cell)
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
const concepts::Sequence< UnderlyingElement > & uelm() const
void recomputeShapefunctions()
const concepts::Z2 mapRho() const
std::string typeOf(const T &t)
Definition output.hh:43
double Real
Definition typedefs.hh:39