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"
16 #include "vectorial/vectorial.hh"
17 #include "hp1D/element.hh"
18 #include "hp2D/quad.hh"
19 
20 namespace 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)
53  {
54  values_ = values;
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 
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:
143 
145 
147  std::unique_ptr<concepts::ShapeFunction1D<Real> > shpfct_;
149  std::unique_ptr<concepts::ShapeFunction1D<Real> > shpfctD_;
160 
161  //Orientation of the mapping, the mapping is not oriented like to topological one in general
163  };
164 
165 } // namespace hp2D
166 
167 #endif // hp2dneumanntrace_hh
A column of a T matrix.
Definition: analytical.hh:18
virtual const concepts::Edge2d & cell() const
Definition: neumannTrace.hh:76
const concepts::Z2 mapRho() const
Returns the orientation of the MappingEdge2d corresponding to the quad.
const concepts::ShapeFunction1D< Real > * shpfct() const
Returns the shape functions.
Definition: neumannTrace.hh:86
A T matrix in sparse notation.
Definition: edgeTest.hh:17
Edge & connector() const
Returns the connector (topology)
Definition: cell1D.hh:59
concepts::Sequence< uint > p_
Polynomial degree along the edge for each side.
virtual ~NTElement_BA()
concepts::Array< Real > values_
The actual storage of the values of the shape functions and their derivatives along the edge.
std::unique_ptr< concepts::ShapeFunction1D< Real > > shpfct_
The shape functions, precomputed on quadrature points.
virtual concepts::Real3d vertex(uint i) const
A 1D cell: edge in 2D.
Definition: cell1D.hh:189
uint n() const
Returns the number of shape functions.
2D hp-FEM for H1-conforming elements.
Real * values_
Values of the shape functions.
Element with cell.
const concepts::Edge2d & cell_
The cell.
uint nP() const
Returns the number of abscissas (in which the shape functions are evaluated)
ShapeFunction1D(const uint n, const uint nP)
Constructor.
virtual const concepts::Edge & support() const
Definition: neumannTrace.hh:73
ushort p() const
concepts::Sequence< uint > dim_
Local number of shape functions coming from one side.
concepts::Sequence< Real > weights_
Weight for each side.
void append(TColumn< F > *T)
Appends the columns to the matrix.
const Real * values() const
Returns the values of the shape functions.
Sequence with operations, output operator, and method of the particular element types.
Definition: sequence.hh:39
virtual const concepts::TMatrix< F > & T() const
Definition: neumannTrace.hh:78
concepts::Z2 mapRho_
concepts::TMatrix< Real > T_
T matrix of the element.
void recomputeShapefunctions()
Recompute shape functions, e.g.
void appendT(concepts::TColumn< F > *T)
Appends the T columns to the T matrix.
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
Definition: neumannTrace.hh:57
const concepts::Sequence< UnderlyingElement > & uelm() const
Returns the Underlying Element(s)
concepts::ElementAndFacette< hp2D::Element< F > > UnderlyingElement
Definition: neumannTrace.hh:42
unsigned short ushort
Abbreviation for unsigned short.
Definition: typedefs.hh:48
Exception class to express a missing feature.
Definition: exceptions.hh:206
Class holding the quadrature rule and the cell of a 1D element.
Definition: element.hh:34
const concepts::ShapeFunction1D< Real > * shpfctD() const
Returns the derivatives of the shape functions.
Definition: neumannTrace.hh:92
ShapeFunction(const uint n, const uint nP, Real *values)
Constructor.
Definition: neumannTrace.hh:51
void addElement(const hp2D::Quad< Real > &quad, uint k, Real weight=1.0)
Adds the contribution to the Neumann trace from the element on one of the (at most two) sides.
Binary group (algebraic): only the values 0 and 1 are represented.
Definition: Zm.hh:16
std::string typeOf(const T &t)
Return the typeid name of a class object.
Definition: output.hh:43
concepts::Array< Real > valuesD_
An edge in the topology.
Definition: topology.hh:73
Element on an edge representing the normal derivative of neighbouring elements, especially their mean...
Definition: neumannTrace.hh:40
std::unique_ptr< concepts::ShapeFunction1D< Real > > shpfctD_
The derivatives of the shape functions.
concepts::Sequence< UnderlyingElement > uelm_
Abstract class for 1D shape function.
NTElement_BA(const concepts::Edge2d &cell)
Constructor.
Container for an element and one facette (edge or face).
Definition: element.hh:113
double Real
Type normally used for a floating point number.
Definition: typedefs.hh:36
virtual std::ostream & info(std::ostream &os) const
Basic namespace for Concepts-2.
Definition: pml_formula.h:16
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich