SOFA API  7dc68564
Open source framework for multi-physics simuation
sofa::geometry::Hexahedron Struct Reference

#include <Hexahedron.h>

Static Public Attributes

static constexpr sofa::Size NumberOfNodes = 8
 
static constexpr ElementType Element_type = ElementType::HEXAHEDRON
 

Public Member Functions

 Hexahedron ()=delete
 

Static Public Member Functions

template<typename Node , typename T = std::decay_t<decltype(*std::begin(std::declval<Node>()))>, typename = std::enable_if_t<std::is_scalar_v<T>>>
static constexpr auto center (const Node &n0, const Node &n1, const Node &n2, const Node &n3, const Node &n4, const Node &n5, const Node &n6, const Node &n7)
 Compute the center of a hexahedron. More...
 
template<typename Node , typename T = std::decay_t<decltype(*std::begin(std::declval<Node>()))>, typename = std::enable_if_t<std::is_scalar_v<T>>>
static constexpr auto barycentricCoefficients (const Node &n0, const Node &n1, const Node &n2, const Node &n3, const Node &n4, const Node &n5, const Node &n6, const Node &n7, const Node &pos)
 Compute the barycentric coefficients of a node in a hexahedron. More...
 
template<typename Node , typename T = std::decay_t<decltype(*std::begin(std::declval<Node>()))>, typename = std::enable_if_t<std::is_scalar_v<T>>>
static constexpr auto squaredDistanceTo (const Node &n0, const Node &n1, const Node &n2, const Node &n3, const Node &n4, const Node &n5, const Node &n6, const Node &n7, const Node &pos)
 Compute the squared distance between a node and the center of a hexahedron. More...
 
template<typename Node , typename T = std::decay_t<decltype(*std::begin(std::declval<Node>()))>, typename = std::enable_if_t<std::is_scalar_v<T>>>
static constexpr auto getPositionFromBarycentricCoefficients (const Node &n0, const Node &n1, const Node &n2, const Node &n3, const Node &n4, const Node &n5, const Node &n6, const Node &n7, const sofa::type::fixed_array< SReal, 3 > &baryC)
 Compute a position from a given set of barycentric coefficients and the associated hexahedron. More...
 
template<typename Node , typename T = std::decay_t<decltype(*std::begin(std::declval<Node>()))>, typename = std::enable_if_t<std::is_scalar_v<T>>>
static constexpr auto volume (const Node &n0, const Node &n1, const Node &n2, const Node &n3, const Node &n4, const Node &n5, const Node &n6, const Node &n7)
 Compute the volume of a hexahedron. More...
 

Attribute details

◆ Element_type

constexpr ElementType sofa::geometry::Hexahedron::Element_type = ElementType::HEXAHEDRON
staticconstexpr

◆ NumberOfNodes

constexpr sofa::Size sofa::geometry::Hexahedron::NumberOfNodes = 8
staticconstexpr

Constructor details

◆ Hexahedron()

sofa::geometry::Hexahedron::Hexahedron ( )
delete

Function details

◆ barycentricCoefficients()

template<typename Node , typename T = std::decay_t<decltype(*std::begin(std::declval<Node>()))>, typename = std::enable_if_t<std::is_scalar_v<T>>>
static constexpr auto sofa::geometry::Hexahedron::barycentricCoefficients ( const Node &  n0,
const Node &  n1,
const Node &  n2,
const Node &  n3,
const Node &  n4,
const Node &  n5,
const Node &  n6,
const Node &  n7,
const Node &  pos 
)
inlinestaticconstexpr

Compute the barycentric coefficients of a node in a hexahedron.

Remarks
Due to some optimizations, the order of nodes given as parameter is necessary.
Template Parameters
Nodeiterable container, with operator[]
Tscalar
Parameters
n0,n1,n2,n3,n4,n5,n6,n7nodes of the hexahedron
posposition of the node from which the coefficients will be computed
Returns
A Vec3 container with the barycentric coefficients of the given node

◆ center()

template<typename Node , typename T = std::decay_t<decltype(*std::begin(std::declval<Node>()))>, typename = std::enable_if_t<std::is_scalar_v<T>>>
static constexpr auto sofa::geometry::Hexahedron::center ( const Node &  n0,
const Node &  n1,
const Node &  n2,
const Node &  n3,
const Node &  n4,
const Node &  n5,
const Node &  n6,
const Node &  n7 
)
inlinestaticconstexpr

Compute the center of a hexahedron.

Remarks
The order of nodes given as parameter is not necessary.
Template Parameters
Nodeiterable container, with operator[]
Tscalar
Parameters
n0,n1,n2,n3,n4,n5,n6,n7nodes of the hexahedron
Returns
Center of the hexahedron (same type as the given nodes)

◆ getPositionFromBarycentricCoefficients()

template<typename Node , typename T = std::decay_t<decltype(*std::begin(std::declval<Node>()))>, typename = std::enable_if_t<std::is_scalar_v<T>>>
static constexpr auto sofa::geometry::Hexahedron::getPositionFromBarycentricCoefficients ( const Node &  n0,
const Node &  n1,
const Node &  n2,
const Node &  n3,
const Node &  n4,
const Node &  n5,
const Node &  n6,
const Node &  n7,
const sofa::type::fixed_array< SReal, 3 > &  baryC 
)
inlinestaticconstexpr

Compute a position from a given set of barycentric coefficients and the associated hexahedron.

Remarks
The order of nodes given as parameter is necessary.
Template Parameters
Nodeiterable container, with operator* applicable with a scalar
Tscalar
Parameters
n0,n1,n2,n3,n4,n5,n6,n7nodes of the hexahedron
baryCbarycentric coefficients
Returns
Position computed from the coefficients, as a Node type

◆ squaredDistanceTo()

template<typename Node , typename T = std::decay_t<decltype(*std::begin(std::declval<Node>()))>, typename = std::enable_if_t<std::is_scalar_v<T>>>
static constexpr auto sofa::geometry::Hexahedron::squaredDistanceTo ( const Node &  n0,
const Node &  n1,
const Node &  n2,
const Node &  n3,
const Node &  n4,
const Node &  n5,
const Node &  n6,
const Node &  n7,
const Node &  pos 
)
inlinestaticconstexpr

Compute the squared distance between a node and the center of a hexahedron.

Remarks
Due to some optimizations, the order of nodes given as parameter is necessary.
Template Parameters
Nodeiterable container, with operator[]
Tscalar
Parameters
n0,n1,n2,n3,n4,n5,n6,n7nodes of the hexahedron
posposition of the node from which the distance will be computed
Returns
Distance from the node and the center of the hexahedron, as a T scalar

◆ volume()

template<typename Node , typename T = std::decay_t<decltype(*std::begin(std::declval<Node>()))>, typename = std::enable_if_t<std::is_scalar_v<T>>>
static constexpr auto sofa::geometry::Hexahedron::volume ( const Node &  n0,
const Node &  n1,
const Node &  n2,
const Node &  n3,
const Node &  n4,
const Node &  n5,
const Node &  n6,
const Node &  n7 
)
inlinestaticconstexpr

Compute the volume of a hexahedron.

Remarks
non optimized version: just return the sum of the 6 inner-tetrahedra
Template Parameters
Nodeiterable container
Tscalar
Parameters
n0,n1,n2,n3,n4,n5,n6,n7nodes of the hexahedron
Returns
Volume of the hexahedron, as a T scalar