hexahedron.hh

Go to the documentation of this file.
1 
6 #ifndef hexahedron_hh
7 #define hexahedron_hh
8 
9 #include <iostream>
10 #include <memory>
11 #include "basics/typedefs.hh"
13 #include "toolbox/array.hh"
14 #include "geometry/cell3D.hh"
15 #include "geometry/integral.hh"
16 #include "integration/quadRule.hh"
18 #include "space/tmatrix.hh"
19 #include "space/hpMethod.hh"
20 #include "element.hh"
21 #include "hexahedronGraphics.hh"
22 
23 #include <map>
24 
25 
26 namespace hp3D {
27  //forward declaration
28  class HexahedronGraphics;
29  using concepts::Real;
30 
31  // ************************************************************ Hexahedron **
32 
37  class Hexahedron : public Element<Real>, public concepts::IntegrationCell {
38  public:
48  virtual ~Hexahedron();
49 
50  virtual const concepts::Hexahedron& support() const {
51  return cell_.connector();
52  }
53  virtual concepts::Real3d vertex(uint i) const {
54  return cell_.vertex(i);
55  }
56  virtual const concepts::Hexahedron3d& cell() const {
57  return cell_;
58  }
59 
61  inline const concepts::Karniadakis<1, 0>* shpfctX() const {
62  return shpfctX_.get();
63  }
65  inline const concepts::Karniadakis<1, 0>* shpfctY() const {
66  return shpfctY_.get();
67  }
69  inline const concepts::Karniadakis<1, 0>* shpfctZ() const {
70  return shpfctZ_.get();
71  }
72 
74  inline const concepts::Karniadakis<1, 1>* shpfctDX() const {
75  return shpfctDX_.get();
76  }
78  inline const concepts::Karniadakis<1, 1>* shpfctDY() const {
79  return shpfctDY_.get();
80  }
82  inline const concepts::Karniadakis<1, 1>* shpfctDZ() const {
83  return shpfctDZ_.get();
84  }
85 
87  inline const concepts::QuadratureRule1d* integrationX() const {
88  return intX_.get();
89  }
91  inline const concepts::QuadratureRule1d* integrationY() const {
92  return intY_.get();
93  }
95  inline const concepts::QuadratureRule1d* integrationZ() const {
96  return intZ_.get();
97  }
98 
101  inline concepts::Real3d chi(const Real x, const Real y,
102  const Real z) const {
103  return cell_.chi(x, y, z);
104  }
105 
107  inline concepts::MapReal3d jacobian(const Real x, const Real y,
108  const Real z) const {
109  concepts::MapReal3d jac = cell_.jacobian(x, y, z);
110  jac *= 0.5;
111  return jac;
112  }
113 
115  inline concepts::MapReal3d jacobianInverse(const Real x, const Real y,
116  const Real z) const {
117  return jacobian(x, y, z).inverse();
118  }
119 
121  inline Real jacobianDeterminant(const Real x, const Real y,
122  const Real z) const {
123  return jacobian(x, y, z).determinant();
124  }
126  inline concepts::MapReal3d hessian(const uint i, const Real x, const Real y,
127  const Real z) const {
128  concepts::MapReal3d hes = cell_.hessian(i,x, y, z);
129  hes *= 0.25;
130  return hes;
131  }
132 
140  void setStrategy(const concepts::Hex3dSubdivision* strategy = 0)
141  throw (concepts::StrategyChange) {
142  cell_.setStrategy(strategy);
143  }
144 
149  return cell_.getStrategy();
150  }
151 
153 
154  virtual bool operator<(const Element<Real>& elm) const;
155 
161  bool hasSameMatrix(const Hexahedron& elm) const;
162 
164  void edgeP(const uint i, const concepts::AdaptiveControlP<1>& p) {
166  edges_[i] = &p;
167  }
169  const concepts::AdaptiveControlP<1>& edgeP(const uint i) const {
172  return *(edges_[i]);
173  }
175  void faceP(const uint i, const concepts::AdaptiveControlP<2>& p) {
177  faces_[i] = &p;
178  }
180  const concepts::AdaptiveControlP<2>& faceP(const uint i) const {
183  return *(faces_[i]);
184  }
185 
194  return rule_;
195  }
196 
197  virtual bool quadraturePoint(uint i, intPoint &p, intFormType form,
198  bool localCoord) const;
199 
201  void recomputeShapefunctions(const uint nq[2]);
202  protected:
203  virtual std::ostream& info(std::ostream& os) const;
206  const concepts::QuadratureRule1d *intY,
207  const concepts::QuadratureRule1d *intZ);
208  private:
212  std::unique_ptr<concepts::Karniadakis<1, 0> > shpfctX_, shpfctY_, shpfctZ_;
214  std::unique_ptr<concepts::Karniadakis<1, 1> > shpfctDX_, shpfctDY_, shpfctDZ_;
216  std::unique_ptr<concepts::QuadratureRule1d> intX_, intY_, intZ_;
217 
222 
223  static std::unique_ptr<HexahedronGraphics> graphics_;
224 
225  // Factory for the integration rule
227 
229  ushort p_[3];
230 
231  };
232 
233 
234  // ****************************************************** ArrayHexaWeights **
235 
238  class ArrayHexaWeights : public concepts::Array<Real> {
239  public:
245  };
246 
247 } // namespace hp3D
248 
249 #endif // hexahedron_hh
const concepts::Karniadakis< 1, 1 > * shpfctDX() const
Returns the derivatives of the shape functions in x direction.
Definition: hexahedron.hh:74
concepts::MapReal3d jacobian(const Real x, const Real y, const Real z) const
Computes the Jacobian.
Definition: hexahedron.hh:107
A column of a T matrix.
Definition: analytical.hh:18
Mapping< F, DimX, DimY > inverse() const
Returns the inverse of the matrix.
const concepts::Karniadakis< 1, 0 > * shpfctY() const
Returns the shape functions in y direction.
Definition: hexahedron.hh:65
std::unique_ptr< concepts::Karniadakis< 1, 0 > > shpfctX_
The shape functions.
Definition: hexahedron.hh:212
void setStrategy(const concepts::Hex3dSubdivision *strategy=0)
Sets the subdivision strategy of the underlying cell of this element.
Definition: hexahedron.hh:140
A 3D FEM element: a hexahedron.
Definition: hexahedron.hh:37
F determinant() const
Returns the determinant of the matrix (only valid for square matrices)
Integration point consisting of coordinates and intermediate data.
Definition: integral.hh:34
const concepts::Karniadakis< 1, 0 > * shpfctZ() const
Returns the shape functions in z direction.
Definition: hexahedron.hh:69
void edgeP(const uint i, const concepts::AdaptiveControlP< 1 > &p)
Set polynomial degree of edge i to p.
Definition: hexahedron.hh:164
ArrayHexaWeights(const Hexahedron &hex)
Constructor for a Hexahedron.
virtual ~Hexahedron()
const concepts::Hex3dSubdivision * getStrategy() const
Returns the subdivision strategy of the underlying cell of this element.
Definition: hexahedron.hh:148
virtual const concepts::Hexahedron & support() const
Definition: hexahedron.hh:50
const concepts::QuadratureRule1d * integrationX() const
Returns the integration rule in x direction.
Definition: hexahedron.hh:87
A 3D cell: hexahedron.
Definition: cell3D.hh:317
const concepts::Karniadakis< 1, 1 > * shpfctDZ() const
Returns the shape functions in z direction.
Definition: hexahedron.hh:82
intFormType
Integration form, which determines terms coming from integration over reference element.
Definition: integral.hh:29
Class for creation of a quadrature rule.
Definition: quadRule.hh:224
std::unique_ptr< concepts::QuadratureRule1d > intZ_
Definition: hexahedron.hh:216
Real jacobianDeterminant(const Real x, const Real y, const Real z) const
Computes the determinant of the Jacobian.
Definition: hexahedron.hh:121
std::unique_ptr< concepts::Karniadakis< 1, 1 > > shpfctDY_
Definition: hexahedron.hh:214
virtual const concepts::ElementGraphics< Real > * graphics() const
void computeShapefunctions_(const concepts::QuadratureRule1d *intX, const concepts::QuadratureRule1d *intY, const concepts::QuadratureRule1d *intZ)
gets the shapefunctions, used in both constructors
virtual bool operator<(const Element< Real > &elm) const
Comparison operator for elements.
#define conceptsAssert(cond, exc)
Assert that a certain condition is fulfilled.
Definition: exceptions.hh:394
void recomputeShapefunctions(const uint nq[2])
static concepts::QuadRuleFactory & rule()
Access to the quadrature rule, which is valid for all elements of this type (hp3D::Hexaedron).
Definition: hexahedron.hh:193
const Hex3dSubdivision * getStrategy() const
Returns the subdivision strategy of this hexahedron.
Definition: cell3D.hh:486
void recomputeShapefunctions()
virtual const concepts::Hexahedron3d & cell() const
Definition: hexahedron.hh:56
Class to represent the quadrature weights on all quadrature points.
Definition: hexahedron.hh:238
const concepts::AdaptiveControlP< 1 > * edges_[12]
Polynomial degree of edges.
Definition: hexahedron.hh:219
concepts::Hexahedron3d & cell_
The cell.
Definition: hexahedron.hh:210
void faceP(const uint i, const concepts::AdaptiveControlP< 2 > &p)
Set polynomial degree of face i to p.
Definition: hexahedron.hh:175
std::unique_ptr< concepts::Karniadakis< 1, 1 > > shpfctDZ_
Definition: hexahedron.hh:214
virtual std::ostream & info(std::ostream &os) const
const concepts::QuadratureRule1d * integrationY() const
Returns the integration rule in y direction.
Definition: hexahedron.hh:91
const ushort * p() const
Returns the polynomial degree.
Definition: element.hh:51
ushort p_[3]
Polynomial degree.
Definition: hexahedron.hh:229
Hexahedron(concepts::Hexahedron3d &cell, const ushort *p, concepts::TColumn< Real > *T0, concepts::TColumn< Real > *T1)
Constructor.
concepts::MapReal3d hessian(const uint i, const Real x, const Real y, const Real z) const
Computes the Hessian.
Definition: hexahedron.hh:126
Real3d vertex(uint i) const
Returns the coordinates of the ith vertex.
static std::unique_ptr< HexahedronGraphics > graphics_
Definition: hexahedron.hh:223
Real3d chi(Real xi, Real eta, Real zeta) const
The element map.
An array of objects.
Definition: bilinearForm.hh:23
Exception class for assertions.
Definition: exceptions.hh:258
A hexahedron in the topology.
Definition: topology3D.hh:134
const concepts::Karniadakis< 1, 0 > * shpfctX() const
Returns the shape functions in x direction.
Definition: hexahedron.hh:61
const concepts::QuadratureRule1d * integrationZ() const
Returns the integration rule in z direction.
Definition: hexahedron.hh:95
virtual concepts::Real3d vertex(uint i) const
Returns the coordinates of the ith vertex of this element.
Definition: hexahedron.hh:53
std::unique_ptr< concepts::Karniadakis< 1, 1 > > shpfctDX_
The derivatives of the shape functions.
Definition: hexahedron.hh:214
const concepts::AdaptiveControlP< 2 > * faces_[6]
Polynomial degree of faces.
Definition: hexahedron.hh:221
concepts::Real3d chi(const Real x, const Real y, const Real z) const
Computes the element map.
Definition: hexahedron.hh:101
MapReal3d jacobian(const Real xi, const Real eta, const Real zeta) const
Computes the Jacobian for xi, eta, zeta .
MapReal3d hessian(const uint i, const Real xi, const Real eta, const Real zeta) const
const concepts::Karniadakis< 1, 1 > * shpfctDY() const
Returns the shape functions in y direction.
Definition: hexahedron.hh:78
virtual bool quadraturePoint(uint i, intPoint &p, intFormType form, bool localCoord) const
Delivers a quadrature point.
Cell over which can be integrated.
Definition: integral.hh:24
std::unique_ptr< concepts::Karniadakis< 1, 0 > > shpfctZ_
Definition: hexahedron.hh:212
const concepts::AdaptiveControlP< 1 > & edgeP(const uint i) const
Get polynomial degree of edge i.
Definition: hexahedron.hh:169
unsigned short ushort
Abbreviation for unsigned short.
Definition: typedefs.hh:48
concepts::MapReal3d jacobianInverse(const Real x, const Real y, const Real z) const
Computes the inverse of the Jacobian.
Definition: hexahedron.hh:115
Abstract class for a 3D FEM element.
Definition: element.hh:29
Hexahedron & connector() const
Returns the connector.
Definition: cell3D.hh:382
Quadrature rule for numerical integration.
Definition: quadRule.hh:30
std::unique_ptr< concepts::Karniadakis< 1, 0 > > shpfctY_
Definition: hexahedron.hh:212
3D hp-FEM for H1-conforming elements.
Definition: meshDX.hh:23
static concepts::QuadRuleFactory rule_
Definition: hexahedron.hh:226
std::unique_ptr< concepts::QuadratureRule1d > intY_
Definition: hexahedron.hh:216
bool hasSameMatrix(const Hexahedron &elm) const
Returns true if element matrix is the same.
Interface for geometrical subdivision strategies for hexahedrons.
Definition: cell3D.hh:159
double Real
Type normally used for a floating point number.
Definition: typedefs.hh:36
void setStrategy(const Hex3dSubdivision *strategy=0)
Sets the subdivision strategy of this hexahedron.
Basic namespace for Concepts-2.
Definition: pml_formula.h:16
const concepts::AdaptiveControlP< 2 > & faceP(const uint i) const
Get polynomial degree of face i.
Definition: hexahedron.hh:180
std::unique_ptr< concepts::QuadratureRule1d > intX_
The integration rules.
Definition: hexahedron.hh:216
Part of the multidimensional expansion bases for the shape functions of Karniadakis and Sherwin.
Definition: karniadakis.hh:163
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich