neumannTraceElement3d.hh

Go to the documentation of this file.
1 
7 #ifndef hp3dntelement3d_hh
8 #define hp3dntelement3d_hh
9 
10 #include "toolbox/sequence.hh"
11 #include "toolbox/sharedPointer.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 
20 namespace hp3D{
21 
22  using concepts::Real;
23 
24  // forward declaration
25  class Hexahedron;
26 
27 
28  // **************************************************** NeumannTraceElement **
29 
41  template<class F = Real>
42  class NeumannTraceElement3d : public hp2D::Quad<F>{
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 
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
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
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 
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
189 
191 
200 
203 
204 
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 
concepts::Array< Real > values_
The actual storage of the values of the shape functions along the edge on integration points.
A 3D FEM element: a hexahedron.
Definition: hexahedron.hh:37
NeumannTraceElement3d(concepts::QuadNd &cell, const ushort *p)
Constructor.
A T matrix in sparse notation.
Definition: edgeTest.hh:17
mapPart(const concepts::Cell2 *coarseCellpart, Real t0, Real t1)
concepts::ElementAndFacette< hp3D::Element< F > > UnderlyingElement
void addElement(const hp3D::Hexahedron &hex, uint k, Real weight=1.0)
Adds the contribution to the Neumann trace from the element on one of the (at most two) sides.
2D hp-FEM for H1-conforming elements.
const concepts::Sequence< UnderlyingElement > & uelm() const
Returns the Underlying Element(s)
void compute_(const concepts::Real2d &q, concepts::Array< Real > &val) const
Routine calculating the normal derivative of the basis functions (of underlying quads),...
Base class for a quadrilateral in any dimension.
Definition: cell2D.hh:187
void buildT()
After adding the Elements this routine builds the TMatrix and polynomial informations on the element.
Sequence with operations, output operator, and method of the particular element types.
Definition: sequence.hh:39
Two dimensional cell.
Definition: cell.hh:89
concepts::Sequence< UnderlyingElement > uelm_
unsigned short ushort
Abbreviation for unsigned short.
Definition: typedefs.hh:48
A 2D FEM element: a quad.
Definition: bf_advection.hh:44
virtual std::ostream & info(std::ostream &os) const
concepts::TMatrix< Real > T_
std::unique_ptr< concepts::ShapeFunction1D< Real > > shpfct_
3D hp-FEM for H1-conforming elements.
Definition: meshDX.hh:23
concepts::Sequence< Real > weights_
Weight for each side. This defines the TraceType on the Element.
void computeShpfctVals(const concepts::Real2d &xP, concepts::Array< Real > &val) const
Evaluation method to get basisfunctions evaluated on requested points this should be used for functio...
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 const concepts::TMatrix< F > & T() const
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich