A 3D FEM element: a hexahedron. More...
#include <hexahedron.hh>


Public Types | |
| typedef F | type |
Public Member Functions | |
| void | appendT (concepts::TColumn< Real > *T) |
| Appends the T columns to the T matrix. | |
| virtual const concepts::Hexahedron3d & | cell () const |
| Returns the cell on which the element is built. | |
| concepts::Real3d | chi (const Real x, const Real y, const Real z) const |
| Computes the element map. | |
| void | edgeP (const uint i, const concepts::AdaptiveControlP< 1 > &p) |
Set polynomial degree of edge i to p. | |
| const concepts::AdaptiveControlP< 1 > & | edgeP (const uint i) const |
Get polynomial degree of edge i. | |
| Real3d | elemMap (const Real coord_local) const |
| Real3d | elemMap (const Real2d &coord_local) const |
| Real3d | elemMap (const Real3d &coord_local) const |
| void | faceP (const uint i, const concepts::AdaptiveControlP< 2 > &p) |
Set polynomial degree of face i to p. | |
| const concepts::AdaptiveControlP< 2 > & | faceP (const uint i) const |
Get polynomial degree of face i. | |
| const concepts::Hex3dSubdivision * | getStrategy () const |
| Returns the subdivision strategy of the underlying cell of this element. | |
| virtual const concepts::ElementGraphics < Real > * | graphics () const |
| bool | hasSameMatrix (const Hexahedron &elm) const |
| Returns true if element matrix is the same. | |
| Hexahedron (concepts::Hexahedron3d &cell, const ushort *p, concepts::TColumn< Real > *T0, concepts::TColumn< Real > *T1) | |
| Constructor. | |
| const concepts::Quadrature < INTEGRATION_RULE > * | integrationX () const |
| Returns the integration rule in x direction. | |
| const concepts::Quadrature < INTEGRATION_RULE > * | integrationY () const |
| Returns the integration rule in y direction. | |
| const concepts::Quadrature < INTEGRATION_RULE > * | integrationZ () const |
| Returns the integration rule in z direction. | |
| concepts::MapReal3d | jacobian (const Real x, const Real y, const Real z) const |
| Computes the Jacobian. | |
| Real | jacobianDeterminant (const Real x, const Real y, const Real z) const |
| Computes the determinant of the Jacobian. | |
| concepts::MapReal3d | jacobianInverse (const Real x, const Real y, const Real z) const |
| Computes the inverse of the Jacobian. | |
| virtual bool | operator< (const Element< F > &elm) const =0 |
| Comparison operator for elements. | |
| virtual bool | operator< (const Element< Real > &elm) const |
| const ushort * | p () const |
| Returns the polynomial degree. | |
| void | setStrategy (const concepts::Hex3dSubdivision *strategy=0) throw (concepts::StrategyChange) |
| Sets the subdivision strategy of the underlying cell of this element. | |
| const concepts::Karniadakis< 1, 1 > * | shpfctDX () const |
| Returns the derivatives of the shape functions in x direction. | |
| const concepts::Karniadakis< 1, 1 > * | shpfctDY () const |
| Returns the shape functions in y direction. | |
| const concepts::Karniadakis< 1, 1 > * | shpfctDZ () const |
| Returns the shape functions in z direction. | |
| const concepts::Karniadakis< 1, 0 > * | shpfctX () const |
| Returns the shape functions in x direction. | |
| const concepts::Karniadakis< 1, 0 > * | shpfctY () const |
| Returns the shape functions in y direction. | |
| const concepts::Karniadakis< 1, 0 > * | shpfctZ () const |
| Returns the shape functions in z direction. | |
| virtual const concepts::Hexahedron & | support () const |
| Returns the topolgical support of the element. | |
| virtual const concepts::TMatrixBase< Real > & | T () const |
| Returns the T matrix of the element. | |
| uint & | tag () |
| Returns the tag. | |
| virtual concepts::Real3d | vertex (uint i) const |
| Returns the coordinates of the ith vertex of this element. | |
| virtual | ~Hexahedron () |
Protected Member Functions | |
| virtual std::ostream & | info (std::ostream &os) const |
| Returns information in an output stream. | |
Protected Attributes | |
| concepts::TMatrix< Real > | T_ |
| T matrix of the element. | |
Private Attributes | |
| concepts::Hexahedron3d & | cell_ |
| The cell. | |
| const concepts::AdaptiveControlP< 1 > * | edges_ [12] |
| Polynomial degree of edges. | |
| const concepts::AdaptiveControlP< 2 > * | faces_ [6] |
| Polynomial degree of faces. | |
| std::auto_ptr < concepts::Quadrature < INTEGRATION_RULE > > | intX_ |
| The integration rules. | |
| std::auto_ptr < concepts::Quadrature < INTEGRATION_RULE > > | intY_ |
| std::auto_ptr < concepts::Quadrature < INTEGRATION_RULE > > | intZ_ |
| std::auto_ptr < concepts::Karniadakis< 1, 1 > > | shpfctDX_ |
| The derivatives of the shape functions. | |
| std::auto_ptr < concepts::Karniadakis< 1, 1 > > | shpfctDY_ |
| std::auto_ptr < concepts::Karniadakis< 1, 1 > > | shpfctDZ_ |
| std::auto_ptr < concepts::Karniadakis< 1, 0 > > | shpfctX_ |
| The shape functions. | |
| std::auto_ptr < concepts::Karniadakis< 1, 0 > > | shpfctY_ |
| std::auto_ptr < concepts::Karniadakis< 1, 0 > > | shpfctZ_ |
Static Private Attributes | |
| static std::auto_ptr < HexahedronGraphics > | graphics_ |
A 3D FEM element: a hexahedron.
Definition at line 36 of file hexahedron.hh.
typedef F concepts::Element< F >::type [inherited] |
Definition at line 53 of file element.hh.
| hp3D::Hexahedron::Hexahedron | ( | concepts::Hexahedron3d & | cell, |
| const ushort * | p, | ||
| concepts::TColumn< Real > * | T0, | ||
| concepts::TColumn< Real > * | T1 | ||
| ) |
Constructor.
| cell | Cell on which the element is defined |
| p | Polynomial degree (can be anisotropic), this array must have at least 3 entries. |
| T0 | Part of the T matrix |
| T1 | Part of the T matrix |
| virtual hp3D::Hexahedron::~Hexahedron | ( | ) | [virtual] |
| void hp3D::Element< F >::appendT | ( | concepts::TColumn< Real > * | T | ) | [inline, inherited] |
Appends the T columns to the T matrix.
Definition at line 60 of file element.hh.
| virtual const concepts::Hexahedron3d& hp3D::Hexahedron::cell | ( | ) | const [inline, virtual] |
Returns the cell on which the element is built.
Possible are tetrahedrons, hexahedron, prims and pyramids.
Implements hp3D::Element< F >.
Definition at line 52 of file hexahedron.hh.
| concepts::Real3d hp3D::Hexahedron::chi | ( | const Real | x, |
| const Real | y, | ||
| const Real | z | ||
| ) | const [inline] |
Computes the element map.
The reference element is the unit cube.
Definition at line 86 of file hexahedron.hh.

| void hp3D::Hexahedron::edgeP | ( | const uint | i, |
| const concepts::AdaptiveControlP< 1 > & | p | ||
| ) | [inline] |
Set polynomial degree of edge i to p.
Definition at line 135 of file hexahedron.hh.

| const concepts::AdaptiveControlP<1>& hp3D::Hexahedron::edgeP | ( | const uint | i | ) | const [inline] |
Get polynomial degree of edge i.
Definition at line 139 of file hexahedron.hh.
| Real3d concepts::ElementWithCell< F >::elemMap | ( | const Real2d & | coord_local | ) | const [inline, inherited] |
Definition at line 87 of file element.hh.
| Real3d concepts::ElementWithCell< F >::elemMap | ( | const Real3d & | coord_local | ) | const [inline, inherited] |
Definition at line 91 of file element.hh.
| Real3d concepts::ElementWithCell< F >::elemMap | ( | const Real | coord_local | ) | const [inline, inherited] |
Definition at line 83 of file element.hh.
| void hp3D::Hexahedron::faceP | ( | const uint | i, |
| const concepts::AdaptiveControlP< 2 > & | p | ||
| ) | [inline] |
Set polynomial degree of face i to p.
Definition at line 144 of file hexahedron.hh.

| const concepts::AdaptiveControlP<2>& hp3D::Hexahedron::faceP | ( | const uint | i | ) | const [inline] |
Get polynomial degree of face i.
Definition at line 148 of file hexahedron.hh.
| const concepts::Hex3dSubdivision* hp3D::Hexahedron::getStrategy | ( | ) | const [inline] |
Returns the subdivision strategy of the underlying cell of this element.
Definition at line 120 of file hexahedron.hh.

| virtual const concepts::ElementGraphics<Real>* hp3D::Hexahedron::graphics | ( | ) | const [virtual] |
Reimplemented from concepts::Element< F >.
| bool hp3D::Hexahedron::hasSameMatrix | ( | const Hexahedron & | elm | ) | const |
Returns true if element matrix is the same.
This method returns true if the element matrix of this element and elm are the same. To find this out, the length of the edges and their angles and the polynomial degrees are compared.
| virtual std::ostream& hp3D::Hexahedron::info | ( | std::ostream & | os | ) | const [protected, virtual] |
Returns information in an output stream.
Reimplemented from hp3D::Element< F >.
| const concepts::Quadrature<INTEGRATION_RULE>* hp3D::Hexahedron::integrationX | ( | ) | const [inline] |
Returns the integration rule in x direction.
Definition at line 75 of file hexahedron.hh.
| const concepts::Quadrature<INTEGRATION_RULE>* hp3D::Hexahedron::integrationY | ( | ) | const [inline] |
Returns the integration rule in y direction.
Definition at line 78 of file hexahedron.hh.
| const concepts::Quadrature<INTEGRATION_RULE>* hp3D::Hexahedron::integrationZ | ( | ) | const [inline] |
Returns the integration rule in z direction.
Definition at line 81 of file hexahedron.hh.
| concepts::MapReal3d hp3D::Hexahedron::jacobian | ( | const Real | x, |
| const Real | y, | ||
| const Real | z | ||
| ) | const [inline] |
Computes the Jacobian.
Definition at line 91 of file hexahedron.hh.

| Real hp3D::Hexahedron::jacobianDeterminant | ( | const Real | x, |
| const Real | y, | ||
| const Real | z | ||
| ) | const [inline] |
Computes the determinant of the Jacobian.
Definition at line 103 of file hexahedron.hh.

| concepts::MapReal3d hp3D::Hexahedron::jacobianInverse | ( | const Real | x, |
| const Real | y, | ||
| const Real | z | ||
| ) | const [inline] |
Computes the inverse of the Jacobian.
Definition at line 98 of file hexahedron.hh.

| virtual bool hp3D::Element< F >::operator< | ( | const Element< F > & | elm | ) | const [pure virtual, inherited] |
Comparison operator for elements.
| virtual bool hp3D::Hexahedron::operator< | ( | const Element< Real > & | elm | ) | const [virtual] |
| const ushort* hp3D::Element< F >::p | ( | ) | const [inline, inherited] |
Returns the polynomial degree.
The returned array has 3 elements.
Definition at line 49 of file element.hh.
| void hp3D::Hexahedron::setStrategy | ( | const concepts::Hex3dSubdivision * | strategy = 0 | ) | throw (concepts::StrategyChange) [inline] |
Sets the subdivision strategy of the underlying cell of this element.
It calls Quad2d::setStrategy.
| strategy | Pointer to an instance of a subdivision strategy. |
| StrategyChange | if the change is not allowed (the change is not allowed if there are children present) |
Definition at line 114 of file hexahedron.hh.

| const concepts::Karniadakis<1,1>* hp3D::Hexahedron::shpfctDX | ( | ) | const [inline] |
Returns the derivatives of the shape functions in x direction.
Definition at line 65 of file hexahedron.hh.
| const concepts::Karniadakis<1,1>* hp3D::Hexahedron::shpfctDY | ( | ) | const [inline] |
Returns the shape functions in y direction.
Definition at line 68 of file hexahedron.hh.
| const concepts::Karniadakis<1,1>* hp3D::Hexahedron::shpfctDZ | ( | ) | const [inline] |
Returns the shape functions in z direction.
Definition at line 71 of file hexahedron.hh.
| const concepts::Karniadakis<1,0>* hp3D::Hexahedron::shpfctX | ( | ) | const [inline] |
Returns the shape functions in x direction.
Definition at line 55 of file hexahedron.hh.
| const concepts::Karniadakis<1,0>* hp3D::Hexahedron::shpfctY | ( | ) | const [inline] |
Returns the shape functions in y direction.
Definition at line 58 of file hexahedron.hh.
| const concepts::Karniadakis<1,0>* hp3D::Hexahedron::shpfctZ | ( | ) | const [inline] |
Returns the shape functions in z direction.
Definition at line 61 of file hexahedron.hh.
| virtual const concepts::Hexahedron& hp3D::Hexahedron::support | ( | ) | const [inline, virtual] |
Returns the topolgical support of the element.
Possible supports for an element are hexahedrons, tetrahedrons, prisms and pyramids.
Implements hp3D::Element< F >.
Definition at line 49 of file hexahedron.hh.

| virtual const concepts::TMatrixBase<Real>& hp3D::Element< F >::T | ( | ) | const [inline, virtual, inherited] |
Returns the T matrix of the element.
Implements concepts::Element< F >.
Definition at line 57 of file element.hh.
| uint& concepts::Element< F >::tag | ( | ) | [inline, inherited] |
Returns the tag.
Definition at line 65 of file element.hh.
| virtual concepts::Real3d hp3D::Hexahedron::vertex | ( | uint | i | ) | const [inline, virtual] |
Returns the coordinates of the ith vertex of this element.
| i | Index of the vertex |
Implements hp3D::Element< F >.
Definition at line 51 of file hexahedron.hh.

concepts::Hexahedron3d& hp3D::Hexahedron::cell_ [private] |
The cell.
Definition at line 156 of file hexahedron.hh.
const concepts::AdaptiveControlP<1>* hp3D::Hexahedron::edges_[12] [private] |
Polynomial degree of edges.
Definition at line 165 of file hexahedron.hh.
const concepts::AdaptiveControlP<2>* hp3D::Hexahedron::faces_[6] [private] |
Polynomial degree of faces.
Definition at line 167 of file hexahedron.hh.
std::auto_ptr<HexahedronGraphics> hp3D::Hexahedron::graphics_ [static, private] |
Definition at line 169 of file hexahedron.hh.
std::auto_ptr<concepts::Quadrature<INTEGRATION_RULE> > hp3D::Hexahedron::intX_ [private] |
The integration rules.
Definition at line 162 of file hexahedron.hh.
std::auto_ptr<concepts::Quadrature<INTEGRATION_RULE> > hp3D::Hexahedron::intY_ [private] |
Definition at line 162 of file hexahedron.hh.
std::auto_ptr<concepts::Quadrature<INTEGRATION_RULE> > hp3D::Hexahedron::intZ_ [private] |
Definition at line 162 of file hexahedron.hh.
std::auto_ptr<concepts::Karniadakis<1,1> > hp3D::Hexahedron::shpfctDX_ [private] |
The derivatives of the shape functions.
Definition at line 160 of file hexahedron.hh.
std::auto_ptr<concepts::Karniadakis<1,1> > hp3D::Hexahedron::shpfctDY_ [private] |
Definition at line 160 of file hexahedron.hh.
std::auto_ptr<concepts::Karniadakis<1,1> > hp3D::Hexahedron::shpfctDZ_ [private] |
Definition at line 160 of file hexahedron.hh.
std::auto_ptr<concepts::Karniadakis<1,0> > hp3D::Hexahedron::shpfctX_ [private] |
The shape functions.
Definition at line 158 of file hexahedron.hh.
std::auto_ptr<concepts::Karniadakis<1,0> > hp3D::Hexahedron::shpfctY_ [private] |
Definition at line 158 of file hexahedron.hh.
std::auto_ptr<concepts::Karniadakis<1,0> > hp3D::Hexahedron::shpfctZ_ [private] |
Definition at line 158 of file hexahedron.hh.
concepts::TMatrix<Real> hp3D::Element< F >::T_ [protected, inherited] |
T matrix of the element.
Definition at line 79 of file element.hh.