neumannTraceElement.hh

Go to the documentation of this file.
1 
8 #ifndef hp2dntelement_hh
9 #define hp2dntelement_hh
10 
11 #include "toolbox/sequence.hh"
12 #include "toolbox/sharedPointer.hh"
13 
15 #include "space/space.hh"
16 #include "space/element.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 // **************************************************** 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  */
56  public:
62  NTShapeFunction(const uint N, const uint nP, Real* values) :
64  values_ = values;
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
100  };
101 
116  const mapPart* coarseCell = 0);
117 
118  //Decontructor
120 
121  virtual const concepts::TMatrix<F>& T() const {
122  return T_;
123  }
124 
133  concepts::Array<Real>& val) const;
134 
138  virtual void recomputeShapefunctions();
139 
143  void buildT();
144 
146  virtual const concepts::ShapeFunction1D<Real>* shpfct() const;
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 
173 protected:
174  virtual std::ostream& info(std::ostream& os) const;
175 
176 private:
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
189 
191 
193 
202 
205 
214  concepts::Array<Real>& val) const;
215 
224  concepts::Array<Real>& val) const;
225 
226 };
227 
228 } // namespace hp2D
229 
230 #endif // hp2dntelement_hh
concepts::Sequence< UnderlyingElement > uelm_
virtual void recomputeShapefunctions()
Recompute shape functions, e.g.
void compute_(const concepts::Array< Real > &q, concepts::Array< Real > &val) const
Routine calculating the normal derivative of the basis functions (of underlying quads),...
A T matrix in sparse notation.
Definition: edgeTest.hh:17
virtual const concepts::EdgeNd & cell() const
Definition: element.hh:99
NeumannTraceElement(const concepts::EdgeNd &cell, uint p, const mapPart *coarseCell=0)
Constructor.
std::unique_ptr< concepts::ShapeFunction1D< Real > > shpfct_
Element on an edge representing the normal derivatives of neighbouring elements, especially their mea...
concepts::Array< Real > values_
The actual storage of the values of the shape functions along the edge on integration points.
uint n() const
Returns the number of shape functions.
virtual const concepts::ShapeFunction1D< Real > * shpfct() const
Returns the shape functions.
2D hp-FEM for H1-conforming elements.
Real * values_
Values of the shape functions.
uint nP() const
Returns the number of abscissas (in which the shape functions are evaluated)
ShapeFunction1D(const uint n, const uint nP)
Constructor.
void computeShpfctVals(const concepts::Array< Real > &xP, concepts::Array< Real > &val) const
Evaluation method to get basisfunctions evaluated on requested points this should be used for functio...
NTShapeFunction(const uint N, const uint nP, Real *values)
Constructor.
mapPart(const concepts::EdgeNd *coarseCellpart, Real t0, Real t1)
virtual const concepts::TMatrix< F > & T() const
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
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.
virtual std::ostream & info(std::ostream &os) const
concepts::TMatrix< Real > T_
void irregularCompute_(const concepts::Array< Real > &q, concepts::Array< Real > &val) const
Routine calculating the normal derivative of the basis functions (of underlying quads),...
A 1D cell in any dimension: edge.
Definition: cell1D.hh:32
void addIrrElement(const hp2D::Quad< Real > &coarseQuad, uint k_coarse, Real weight)
Basically the same as addElement, but with the difference in the meaning.
concepts::ElementAndFacette< hp2D::Element< F > > UnderlyingElement
const concepts::Sequence< UnderlyingElement > & uelm() const
Returns the Underlying Element(s)
std::string typeOf(const T &t)
Return the typeid name of a class object.
Definition: output.hh:43
void buildT()
After adding the Elements this routine builds the TMatrix and polynomial informations on the element.
ushort p() const
Definition: element.hh:103
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
Abstract class for 1D shape function.
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
concepts::Sequence< Real > weights_
Weight for each side. This defines the TraceType on the Element.
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