Class documentation of Concepts

Loading...
Searching...
No Matches
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"
18#include "space/tmatrix.hh"
19#include "space/hpMethod.hh"
20#include "element.hh"
21#include "hexahedronGraphics.hh"
22
23#include <map>
24
25
26namespace 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
88 return intX_.get();
89 }
92 return intY_.get();
93 }
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 {
142 cell_.setStrategy(strategy);
143 }
144
149 return cell_.getStrategy();
150 }
151
152 virtual const concepts::ElementGraphics<Real>* graphics() const;
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
200 void recomputeShapefunctions();
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
219 const concepts::AdaptiveControlP<1>* edges_[12];
221 const concepts::AdaptiveControlP<2>* faces_[6];
222
223 static std::unique_ptr<HexahedronGraphics> graphics_;
224
225 // Factory for the integration rule
226 static concepts::QuadRuleFactory 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
Hexahedron & connector() const
Returns the connector.
Definition cell3D.hh:382
MapReal3d jacobian(const Real xi, const Real eta, const Real zeta) const
Real3d vertex(uint i) const
Returns the coordinates of the ith vertex.
void setStrategy(const Hex3dSubdivision *strategy=0)
Real3d chi(Real xi, Real eta, Real zeta) const
const Hex3dSubdivision * getStrategy() const
Definition cell3D.hh:485
Mapping< F, DimX, DimY > inverse() const
Returns the inverse of the matrix.
F determinant() const
Returns the determinant of the matrix (only valid for square matrices)
ArrayHexaWeights(const Hexahedron &hex)
const ushort * p() const
Definition element.hh:51
const concepts::QuadratureRule1d * integrationZ() const
Returns the integration rule in z direction.
Definition hexahedron.hh:95
virtual concepts::Real3d vertex(uint i) const
Definition hexahedron.hh:53
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
void edgeP(const uint i, const concepts::AdaptiveControlP< 1 > &p)
Set polynomial degree of edge i to p.
const concepts::QuadratureRule1d * integrationX() const
Returns the integration rule in x direction.
Definition hexahedron.hh:87
virtual const concepts::Hexahedron3d & cell() const
Definition hexahedron.hh:56
Hexahedron(concepts::Hexahedron3d &cell, const ushort *p, concepts::TColumn< Real > *T0, concepts::TColumn< Real > *T1)
concepts::MapReal3d jacobian(const Real x, const Real y, const Real z) const
Computes the Jacobian.
void computeShapefunctions_(const concepts::QuadratureRule1d *intX, const concepts::QuadratureRule1d *intY, const concepts::QuadratureRule1d *intZ)
gets the shapefunctions, used in both constructors
virtual bool quadraturePoint(uint i, intPoint &p, intFormType form, bool localCoord) const
const concepts::Karniadakis< 1, 1 > * shpfctDZ() const
Returns the shape functions in z direction.
Definition hexahedron.hh:82
virtual const concepts::Hexahedron & support() const
Definition hexahedron.hh:50
concepts::MapReal3d hessian(const uint i, const Real x, const Real y, const Real z) const
Computes the Hessian.
const concepts::AdaptiveControlP< 2 > & faceP(const uint i) const
Get polynomial degree of face i.
concepts::Real3d chi(const Real x, const Real y, const Real z) const
const concepts::AdaptiveControlP< 1 > & edgeP(const uint i) const
Get polynomial degree of edge i.
const concepts::Karniadakis< 1, 0 > * shpfctY() const
Returns the shape functions in y direction.
Definition hexahedron.hh:65
bool hasSameMatrix(const Hexahedron &elm) const
static concepts::QuadRuleFactory & rule()
Real jacobianDeterminant(const Real x, const Real y, const Real z) const
Computes the determinant of the Jacobian.
const concepts::Karniadakis< 1, 1 > * shpfctDY() const
Returns the shape functions in y direction.
Definition hexahedron.hh:78
const concepts::Karniadakis< 1, 1 > * shpfctDX() const
Returns the derivatives of the shape functions in x direction.
Definition hexahedron.hh:74
const concepts::Karniadakis< 1, 0 > * shpfctZ() const
Returns the shape functions in z direction.
Definition hexahedron.hh:69
const concepts::Karniadakis< 1, 0 > * shpfctX() const
Returns the shape functions in x direction.
Definition hexahedron.hh:61
concepts::MapReal3d jacobianInverse(const Real x, const Real y, const Real z) const
Computes the inverse of the Jacobian.
void setStrategy(const concepts::Hex3dSubdivision *strategy=0)
void faceP(const uint i, const concepts::AdaptiveControlP< 2 > &p)
Set polynomial degree of face i to p.
const concepts::Hex3dSubdivision * getStrategy() const
const concepts::QuadratureRule1d * integrationY() const
Returns the integration rule in y direction.
Definition hexahedron.hh:91
#define conceptsAssert(cond, exc)
double Real
Definition typedefs.hh:39
Definition meshDX.hh:23