Caribou
|
Classes | |
struct | GaussNode |
Public Member Functions | |
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... | |
template<typename MatrixType > | |
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 = traits<Derived>::CanonicalDimension |
static constexpr auto | Dimension = traits<Derived>::Dimension |
static constexpr auto | NumberOfNodesAtCompileTime = traits<Derived>::NumberOfNodesAtCompileTime |
static constexpr auto | NumberOfGaussNodesAtCompileTime = traits<Derived>::NumberOfGaussNodesAtCompileTime |
|
inline |
Construct and return the given boundary element.
Example:
|
inline |
Get the list of node indices of the boundary elements.
The return type is up to the implementation but should normally be:
|
inline |
Return true if the element contains the point located at the given local coordinates.
coordinates | Local coordinates of a point |
eps | If the given point is located barely outside the element, which is, less than this eps value, returns true. |
|
inline |
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.
Example:
|
inline |
Interpolate a value at local coordinates from the given interpolation node values.
The values at nodes must be provided in an Eigen::Matrix where the row i of the matrix is the value (as a row-vector or a scalar) at the node i.
MatrixType | Type of the Eigen matrix containing the values at nodes. |
|
inline |
Compute the Jacobian matrix of the transformation T(xi)-> x evaluated at local coordinates xi.
The Jacobian is defined as:
* * |dx| |du| * |dy| = J |dv| * |dz| |dw| * * 1D canonical element * -------------------- * * 1D manifold: J(u) = dx/du = sum_i dNi/du * x_i * det(J) = |J| * * 2D manifold: J(u) = | dx/du | = | sum dNi/du x_i | * | dy/du | = | sum dNi/du y_i | * det(J) = sqrt(J.dot(J)) * * | dx/du | = | sum dNi/du x_i | * 3D manifold: J(u) = | dy/du | = | sum dNi/du y_i | * | dz/du | = | sum dNi/du z_i | * det(J) = sqrt(J.dot(J)) * * * 2D canonical element * -------------------- * * 1D manifold: N/A * * 2D manifold: J(u,v) = | dx/du dx/dv | | sum dNi/du x_i sum dNi/dv x_i | * | dy/du dy/dv | = | sum dNi/du y_i sum dNi/dv y_i | * det(J) = |det(J)| * * | dx/du dx/dv | | sum dNi/du x_i sum dNi/dv x_i | * 3D manifold: J(u,v) = | dy/du dy/dv | = | sum dNi/du y_i sum dNi/dv y_i | * | dz/du dz/dv | = | sum dNi/du z_i sum dNi/dv z_i | * det(J) = sqrt((J.transpose() * J).determinant()) * = J.col(0).cross(J.col(1)).norm(); * * * 3D canonical element * -------------------- * * 1D manifold: N/A * 2D manifold: N/A * | dx/du dx/dv dx/dw | | sum dNi/du x_i sum dNi/dv x_i sum dNi/dw x_i | * 3D manifold: J(u,v,w) = | dy/du dy/dv dy/dw | = | sum dNi/du y_i sum dNi/dv y_i sum dNi/dw y_i | * | dz/du dz/dv dz/dw | | sum dNi/du z_i sum dNi/dv z_i sum dNi/dw z_i | * * * * * where dNi/du (resp. dv and dw) is the partial derivative of the shape function at node i * w.r.t the local frame of the canonical element evaluated at local coordinate {u, v, w} and * where {xi, yi and zi} are the world coordinates of the position of node i on its element manifold. *
Example:
|
inline |
Get the Lagrange polynomial values evaluated at local coordinates xi w.r.t each element's interpolation nodes.
Example:
|
inline |
Get the local coordinates of a point from its world coordinates by doing a set of Newton-Raphson iterations.
|
inline |
Get the local coordinates of a point from its world coordinates by doing a set of Newton-Raphson iterations.
By taking the Taylor expansion of the transformation with being the world coordinates of a point and its local coordinates within the element, we have
where partial derivatives are evaluated at . We can reformulate with the following matrix form
where , and is the Jacobian of the transformation . Since we are trying to find , we can rearange the last equation to get
Hence, starting from an initial guess at local coordinates , we have the following iterative method:
The iterations stop when
coordinates | The world coordinates of the point from which we want to get the local coordinates. |
starting_point | An approximation of the real local coordinates we want the get. The closer it is to the solution, the faster the Newton will converge. |
residual_tolerance | The threshold of relative norm of the residual at which point the Newton is said to converge (|R|/|R0| < threshold). |
maximum_number_of_iterations | The maximum number of Newton-Raphson iterations we can do before divergence. |
|
inline |
Get the number of boundary elements (ex.
number of triangles in a tetrahedron)