Mesh Oriented datABase
(version 5.4.1)
Array-based unstructured mesh datastructure
|
Class representing a 3-D mapping function (e.g. shape function for volume element) More...
Public Member Functions | |
virtual CartVect | center_xi () const =0 |
Return \(\vec \xi\) corresponding to logical center of element. | |
virtual CartVect | evaluate (const CartVect &xi) const =0 |
Evaluate mapping function (calculate \(\vec x = F($\vec \xi)\) ) | |
virtual Matrix3 | jacobian (const CartVect &xi) const =0 |
Evaluate Jacobian of mapping function. | |
bool | solve_inverse (const CartVect &x, CartVect &xi, double tol) const |
Evaluate inverse of mapping function (calculate \(\vec \xi = F^-1($\vec x)\) ) |
Class representing a 3-D mapping function (e.g. shape function for volume element)
Definition at line 16 of file ElemUtil.cpp.
virtual CartVect moab::ElemUtil::VolMap::center_xi | ( | ) | const [pure virtual] |
Return \(\vec \xi\) corresponding to logical center of element.
Implemented in moab::ElemUtil::LinearHexMap.
Referenced by solve_inverse().
virtual CartVect moab::ElemUtil::VolMap::evaluate | ( | const CartVect & | xi | ) | const [pure virtual] |
Evaluate mapping function (calculate \(\vec x = F($\vec \xi)\) )
Implemented in moab::ElemUtil::LinearHexMap.
Referenced by solve_inverse().
virtual Matrix3 moab::ElemUtil::VolMap::jacobian | ( | const CartVect & | xi | ) | const [pure virtual] |
Evaluate Jacobian of mapping function.
Implemented in moab::ElemUtil::LinearHexMap.
Referenced by solve_inverse().
bool moab::ElemUtil::VolMap::solve_inverse | ( | const CartVect & | x, |
CartVect & | xi, | ||
double | tol | ||
) | const |
Evaluate inverse of mapping function (calculate \(\vec \xi = F^-1($\vec x)\) )
Definition at line 29 of file ElemUtil.cpp.
References center_xi(), moab::Matrix3::determinant(), evaluate(), moab::Matrix3::inverse(), and jacobian().
Referenced by moab::ElemUtil::nat_coords_trilinear_hex().
{ const double error_tol_sqr = tol * tol; double det; xi = center_xi(); CartVect delta = evaluate( xi ) - x; Matrix3 J; while( delta % delta > error_tol_sqr ) { J = jacobian( xi ); det = J.determinant(); if( det < std::numeric_limits< double >::epsilon() ) return false; xi -= J.inverse() * delta; delta = evaluate( xi ) - x; } return true; }