Caribou
Public Types | Public Member Functions | Static Public Attributes | Friends | List of all members
caribou::geometry::Hexahedron< Quadratic > Struct Reference

Detailed Description

Quadratic Hexahedron (20 nodes)

*         v
*         ^
*         |
*  3-----10---2
*  |\     |   |\
*  | 19   |   | 18
* 11  \   |   9  \
*  |   7----14+---6
*  |   |  +-- |-- | ---> u
*  0---+-8-\--1   |
*   \  15   \  \  13
*   16 |     \  17|
*     \|      w  \|
*      4----12----5   u
* 

#include <Hexahedron.h>

Inheritance diagram for caribou::geometry::Hexahedron< Quadratic >:
caribou::geometry::BaseHexahedron< Derived > caribou::geometry::Element< Derived >

Public Types

using Base = BaseHexahedron< Hexahedron< Quadratic > >
 
using LocalCoordinates = typename Base::LocalCoordinates
 
using WorldCoordinates = typename Base::WorldCoordinates
 
using GaussNode = typename Base::GaussNode
 
template<UNSIGNED_INTEGER_TYPE Dim>
using Vector = typename Base::template Vector< Dim >
 
template<UNSIGNED_INTEGER_TYPE Rows, UNSIGNED_INTEGER_TYPE Cols>
using Matrix = typename Base::template Matrix< Rows, Cols >
 
- Public Types inherited from caribou::geometry::BaseHexahedron< Derived >
using Base = Element< Derived >
 
using LocalCoordinates = typename Base::LocalCoordinates
 
using WorldCoordinates = typename Base::WorldCoordinates
 
using GaussNode = typename Base::GaussNode
 
template<UNSIGNED_INTEGER_TYPE Dim>
using Vector = typename Base::template Vector< Dim >
 
template<UNSIGNED_INTEGER_TYPE Rows, UNSIGNED_INTEGER_TYPE Cols>
using Matrix = typename Base::template Matrix< Rows, Cols >
 
- Public Types inherited from caribou::geometry::Element< Derived >
using Scalar = FLOATING_POINT_TYPE
 
using Vector = Eigen::Matrix< Scalar, Dim, 1 >
 
using Matrix = Eigen::Matrix< Scalar, Rows, Cols, Options >
 
using MatrixI = Eigen::Matrix< Scalar, Rows, Cols, Options >
 
using LocalCoordinates = Vector< CanonicalDimension >
 
using WorldCoordinates = Vector< Dimension >
 

Public Member Functions

 Hexahedron (const Hexahedron< Linear > &linear_Hexahedron)
 
 Hexahedron (WorldCoordinates &p0, WorldCoordinates &p1, WorldCoordinates &p2, WorldCoordinates &p3, WorldCoordinates &p4, WorldCoordinates &p5, WorldCoordinates &p6, WorldCoordinates &p7)
 Construct a quadratic Hexahedron from eight nodes.
 
 Hexahedron (const WorldCoordinates &p0, const WorldCoordinates &p1, const WorldCoordinates &p2, const WorldCoordinates &p3, const WorldCoordinates &p4, const WorldCoordinates &p5, const WorldCoordinates &p6, const WorldCoordinates &p7)
 Construct a quadratic Hexahedron from eight nodes.
 
auto edges () const
 Get the list of node indices of the edges. More...
 
- Public Member Functions inherited from caribou::geometry::BaseHexahedron< Derived >
 BaseHexahedron ()=default
 Default empty constructor.
 
template<typename EigenType , REQUIRES(EigenType::RowsAtCompileTime==NumberOfNodesAtCompileTime) >
 BaseHexahedron (Eigen::EigenBase< EigenType > &nodes)
 Constructor from an Eigen matrix containing the positions of the hexahedron's nodes.
 
template<typename EigenType , REQUIRES(EigenType::RowsAtCompileTime==NumberOfNodesAtCompileTime) >
 BaseHexahedron (const Eigen::EigenBase< EigenType > &nodes)
 Constructor from an Eigen matrix containing the positions of the hexahedron's nodes.
 
template<typename ... Nodes>
 BaseHexahedron (const WorldCoordinates &first_node, Nodes &&...remaining_nodes)
 Constructor from a serie of nodes.
 
auto faces () const
 Get the list of node indices of the faces. More...
 
auto frame () const -> Matrix< 3, 3 >
 Extract the orthogonal frame of the element by computing the cross product of the unit vectors from the center position to its projection on opposite faces. More...
 
auto frame (const LocalCoordinates &local_point) const -> Matrix< 3, 3 >
 Extract the frame positioned at the given position (in local coordinates) on the hexa by computing the cross product of the unit vectors from the given position its projection on opposite faces. More...
 
- Public Member Functions inherited from caribou::geometry::Element< Derived >
auto number_of_nodes () const -> UNSIGNED_INTEGER_TYPE
 Get the number of nodes in the element.
 
auto number_of_gauss_nodes () const -> UNSIGNED_INTEGER_TYPE
 Get the number of gauss nodes in the element.
 
auto node (const UNSIGNED_INTEGER_TYPE &index) const
 Get the Node at given index.
 
auto nodes () const
 Get the set of nodes.
 
auto gauss_node (const UNSIGNED_INTEGER_TYPE &index) const -> const GaussNode &
 Get the gauss node at given index.
 
auto gauss_nodes () const -> const std::vector< GaussNode > &
 Get the set of gauss nodes.
 
auto number_of_boundary_elements () const
 Get the number of boundary elements (ex. More...
 
auto boundary_elements_node_indices () const -> const auto &
 Get the list of node indices of the boundary elements. More...
 
auto boundary_element (const UNSIGNED_INTEGER_TYPE &boundary_id) const
 Construct and return the given boundary element. More...
 
auto L (const LocalCoordinates &xi) const -> Vector< NumberOfNodesAtCompileTime >
 Get the Lagrange polynomial values evaluated at local coordinates xi w.r.t each element's interpolation nodes. More...
 
auto dL (const LocalCoordinates &xi) const -> Matrix< NumberOfNodesAtCompileTime, CanonicalDimension >
 Get the Lagrange polynomial derivatives w.r.t the local frame {dL/du} evaluated at local coordinates {u} w.r.t each segment's interpolation nodes. More...
 
auto center () const -> WorldCoordinates
 Get the position at the center of the element.
 
auto world_coordinates (const LocalCoordinates &coordinates) const -> WorldCoordinates
 Get the world coordinates of a point from its local coordinates.
 
auto local_coordinates (const WorldCoordinates &coordinates) const -> LocalCoordinates
 Get the local coordinates of a point from its world coordinates by doing a set of Newton-Raphson iterations. More...
 
auto local_coordinates (const WorldCoordinates &coordinates, const LocalCoordinates &starting_point, const FLOATING_POINT_TYPE &residual_tolerance, const UNSIGNED_INTEGER_TYPE &maximum_number_of_iterations) const -> LocalCoordinates
 Get the local coordinates of a point from its world coordinates by doing a set of Newton-Raphson iterations. More...
 
auto contains_local (const LocalCoordinates &xi, const FLOATING_POINT_TYPE &eps=1e-10) const -> bool
 Return true if the element contains the point located at the given local coordinates. More...
 
auto interpolate (const LocalCoordinates &coordinates, const Eigen::MatrixBase< MatrixType > &values) const
 Interpolate a value at local coordinates from the given interpolation node values. More...
 
auto jacobian (const LocalCoordinates &coordinates) const -> Matrix< Dimension, CanonicalDimension >
 Compute the Jacobian matrix of the transformation T(xi)-> x evaluated at local coordinates xi. More...
 

Static Public Attributes

static constexpr auto CanonicalDimension = Base::CanonicalDimension
 
static constexpr auto Dimension = Base::Dimension
 
static constexpr auto NumberOfNodesAtCompileTime = Base::NumberOfNodesAtCompileTime
 
static constexpr auto NumberOfGaussNodesAtCompileTime = Base::NumberOfGaussNodesAtCompileTime
 
static constexpr FLOATING_POINT_TYPE canonical_nodes [NumberOfNodesAtCompileTime][CanonicalDimension]
 
- Static Public Attributes inherited from caribou::geometry::BaseHexahedron< Derived >
static constexpr auto CanonicalDimension = Base::CanonicalDimension
 
static constexpr auto Dimension = Base::Dimension
 
static constexpr auto NumberOfNodesAtCompileTime = Base::NumberOfNodesAtCompileTime
 
static constexpr auto NumberOfGaussNodesAtCompileTime = Base::NumberOfGaussNodesAtCompileTime
 
- Static Public Attributes inherited from caribou::geometry::Element< Derived >
static constexpr auto CanonicalDimension
 
static constexpr auto Dimension
 
static constexpr auto NumberOfNodesAtCompileTime
 
static constexpr auto NumberOfGaussNodesAtCompileTime
 

Friends

struct Element< Hexahedron< Quadratic > >
 
struct BaseHexahedron< Hexahedron< Quadratic > >
 

Additional Inherited Members

- Protected Attributes inherited from caribou::geometry::BaseHexahedron< Derived >
Matrix< NumberOfNodesAtCompileTime, Dimension > p_nodes
 

Member Function Documentation

◆ edges()

auto caribou::geometry::Hexahedron< Quadratic >::edges ( ) const
inline

Get the list of node indices of the edges.

See also
Element::boundary_elements_node_indices

Member Data Documentation

◆ canonical_nodes

constexpr FLOATING_POINT_TYPE caribou::geometry::Hexahedron< Quadratic >::canonical_nodes[NumberOfNodesAtCompileTime][CanonicalDimension]
staticconstexpr
Initial value:
= {
{-1, -1, -1},
{+1, -1, -1},
{+1, +1, -1},
{-1, +1, -1},
{-1, -1, +1},
{+1, -1, +1},
{+1, +1, +1},
{-1, +1, +1},
{ 0, -1, -1},
{+1, 0, -1},
{ 0, +1, -1},
{-1, 0, -1},
{ 0, -1, +1},
{+1, 0, +1},
{ 0, +1, +1},
{-1, 0, +1},
{-1, -1, 0},
{+1, -1, 0},
{+1, +1, 0},
{-1, +1, 0},
}

The documentation for this struct was generated from the following file: