imaging::FemKernel< fem_types > Class Template Reference
[Finite Element Module]

Computes the values of the shape functions and the element transformation in the integration nodes. More...

#include <FemKernel.hpp>

List of all members.

Public Member Functions

 FemKernel (const Grid< fem_types > &grid)
const Grid< fem_types > & grid () const
void set_element (std::size_t element)
void set_boundary_element (std::size_t boundary_element)
void lazy_set_element (std::size_t element)
float_t shape_value (std::size_t integrator_node, std::size_t element_node) const
const ublas::fixed_vector
< float_t,
fem_types::data_dimension > & 
shape_gradient (std::size_t integrator_node, std::size_t element_node) const
float_t transform_determinant (std::size_t integrator_node) const
float_t shape_boundary_value (std::size_t integrator_node, std::size_t element_node) const
float_t shape_boundary_derivative (std::size_t integrator_node, std::size_t element_node) const
float_t boundary_transform_determinant (std::size_t integrator_node) const
const ublas::fixed_vector
< float_t,
fem_types::data_dimension > & 
boundary_normal () const
std::size_t current_element () const
std::size_t current_boundary_element () const
bool is_boundary_node (std::size_t element_node) const
const ublas::fixed_vector
< float_t,
fem_types::data_dimension > 
boundary_normal_at_node (std::size_t element_node) const


Detailed Description

template<class fem_types>
class imaging::FemKernel< fem_types >

Computes the values of the shape functions and the element transformation in the integration nodes.

FemKernel is a class which is used to evaluate functions in the integration nodes on a given element. The integrator, the shape function and the element type is specified by the template parameter fem_types. If a FemKernel is set to an element by calling set_element() or set_boundary_element(), it precomputes the values and derivatives of the shape functions and the element transformation in the integration points. Then you can use the accessor functions to query these values.

A FemKernel object is automatically created by Assembler::assemble(), Assembler::assemble_stiffness_matrix() or Assembler::assemble_force_vector() and then iteratively set to the elements in a given grid during the assembly process. The initialized kernel is then passed to the equation object of the current FE computation. In the implementation of the equation (cf. EquationInterface) you can query the kernel to compute the equation coefficients which you then pass back to the assembly function.

See also:
EquationInterface

Constructor & Destructor Documentation

template<class fem_types>
imaging::FemKernel< fem_types >::FemKernel ( const Grid< fem_types > &  grid  )  [inline]

Construct kernel from grid. The kernel object stores only a reference to the grid. When the kernel is later on set to an element by calling void set_element(std::size_t element) or set_boundary_element(std::size_t element) the argument element refers to grid.


Member Function Documentation

template<class fem_types>
const Grid<fem_types>& imaging::FemKernel< fem_types >::grid (  )  const [inline]

template<class fem_types>
void imaging::FemKernel< fem_types >::set_element ( std::size_t  element  )  [inline]

Initializes the kernel to element, computes

  • the values and derivatives of the shape functions and
  • the element transform determinant (Jacobian of the transform), and stores them. This function is relatively expensive whereas the retrieval of the computed values later on is cheap. As long as set_element() or lazy_set_element() are not called with a different parameter element, current_element() will return element.

References imaging::inverse().

Referenced by imaging::Assembler::assemble(), imaging::Assembler::assemble_force_vector(), imaging::Assembler::assemble_stiffness_matrix(), imaging::average(), imaging::integrate(), imaging::maximum(), and imaging::minimum().

template<class fem_types>
void imaging::FemKernel< fem_types >::set_boundary_element ( std::size_t  boundary_element  )  [inline]

Initializes the kernel to boundary_element, computes

  • the values and derivatives of the boundary shape functions and
  • the boundary element transform determinant (determinant of the Jacobian of the transform), and stores them. This function is relatively expensive whereas the retrieval of the computed values later on is cheap. As long as set_boundary_element() is not called with a different parameter element, current_boundary_element() will return boundary_element.

References imaging::FemKernel< fem_types >::boundary_normal(), and imaging::inverse().

Referenced by imaging::Assembler::assemble(), imaging::Assembler::assemble_force_vector(), and imaging::Assembler::assemble_stiffness_matrix().

template<class fem_types>
void imaging::FemKernel< fem_types >::lazy_set_element ( std::size_t  element  )  [inline]

Initializes the kernel to element, but does not compute new values of the shape and transformation functions. Use this function if you know for sure that the geometry of the new element is exactly the same as the one of the previous element. As long as set_element() or lazy_set_element() are not called with a different parameter element, current_element() will return element.

Referenced by imaging::Assembler::assemble(), imaging::Assembler::assemble_force_vector(), imaging::Assembler::assemble_stiffness_matrix(), imaging::average(), imaging::integrate(), imaging::maximum(), and imaging::minimum().

template<class fem_types>
float_t imaging::FemKernel< fem_types >::shape_value ( std::size_t  integrator_node,
std::size_t  element_node 
) const [inline]

Returns the value of the shape function assciated with element_node in the integration node integrator_node.

Referenced by imaging::SimpleEquationAdaptor< equation_t >::force_vector(), imaging::Grid< fem_types >::interpolate_value(), and imaging::SimpleEquationAdaptor< equation_t >::stiffness_matrix().

template<class fem_types>
const ublas::fixed_vector<float_t, fem_types::data_dimension>& imaging::FemKernel< fem_types >::shape_gradient ( std::size_t  integrator_node,
std::size_t  element_node 
) const [inline]

Returns the gradient of the shape function assciated with element_node in the integration node integrator_node.

Referenced by imaging::SimpleEquationAdaptor< equation_t >::force_vector(), imaging::Grid< fem_types >::interpolate_gradient(), and imaging::SimpleEquationAdaptor< equation_t >::stiffness_matrix().

template<class fem_types>
float_t imaging::FemKernel< fem_types >::transform_determinant ( std::size_t  integrator_node  )  const [inline]

Returns the value of the transformation determinant (determinant of the Jacobian of the transform) in the integration node integrator_node. This value reflects the distortion of the actual element with respect to its reference element. Note that integrator_node refers to the corresponding node of the boundary integrator.

Referenced by imaging::Assembler::assemble(), imaging::Assembler::assemble_force_vector(), imaging::Assembler::assemble_stiffness_matrix(), imaging::average(), and imaging::integrate().

template<class fem_types>
float_t imaging::FemKernel< fem_types >::shape_boundary_value ( std::size_t  integrator_node,
std::size_t  element_node 
) const [inline]

Returns the value of shape function associated with element_node on the parent element of the current boundary element in the integration node integrator_node on the boundary element. This means that the shape function on the element is evaluated at the boundary of the element. Keep in mind that here the element is not the current element but the parent element of the current boundary element. I.e. this function depends only on the current boundary element. Note that integrator_node refers to the corresponding node of the boundary integrator.

Referenced by imaging::SimpleEquationAdaptor< equation_t >::force_vector_at_boundary(), imaging::Grid< fem_types >::interpolate_boundary_value(), and imaging::SimpleEquationAdaptor< equation_t >::stiffness_matrix_at_boundary().

template<class fem_types>
float_t imaging::FemKernel< fem_types >::shape_boundary_derivative ( std::size_t  integrator_node,
std::size_t  element_node 
) const [inline]

Returns the gradient of shape function associated with element_node on the parent element of the current boundary element in the integration node integrator_node on the boundary element. This means that the gradient of the shape function on the element is evaluated at the boundary of the element. Keep in mind that here the element is not the current element but the parent element of the current boundary element. I.e. this function depends only on the current boundary element. Note that integrator_node refers to the corresponding node of the boundary integrator.

Referenced by imaging::Grid< fem_types >::interpolate_boundary_derivative(), and imaging::SimpleEquationAdaptor< equation_t >::stiffness_matrix_at_boundary().

template<class fem_types>
float_t imaging::FemKernel< fem_types >::boundary_transform_determinant ( std::size_t  integrator_node  )  const [inline]

Returns the value of the transformation determinant (determinant of the Jacobian of the transform) of the transformation which transforms the current boundary element to the reference element. This value reflects the distortion of the actual boundary element with respect to its reference element.

Referenced by imaging::Assembler::assemble(), imaging::Assembler::assemble_force_vector(), imaging::Assembler::assemble_stiffness_matrix(), and imaging::SimpleEquationAdaptor< equation_t >::stiffness_matrix_at_boundary().

template<class fem_types>
const ublas::fixed_vector<float_t, fem_types::data_dimension>& imaging::FemKernel< fem_types >::boundary_normal (  )  const [inline]

Returns the boundary normal in the integration node integrator_node of the current boundary element. Note that integrator_node refers to the corresponding node of the boundary integrator.

Referenced by imaging::FemKernel< fem_types >::set_boundary_element().

template<class fem_types>
std::size_t imaging::FemKernel< fem_types >::current_element (  )  const [inline]

Returns the index of the current element. In most situations you should assume that the kernel is initialized to this element.

Referenced by imaging::Grid< fem_types >::global_coordinates(), imaging::Grid< fem_types >::interpolate_gradient(), and imaging::Grid< fem_types >::interpolate_value().

template<class fem_types>
std::size_t imaging::FemKernel< fem_types >::current_boundary_element (  )  const [inline]

Returns the index of the current boundary element. In most situations you should assume that the kernel is initialized to this boundary element.

Referenced by imaging::Grid< fem_types >::interpolate_boundary_derivative(), and imaging::Grid< fem_types >::interpolate_boundary_value().

template<class fem_types>
bool imaging::FemKernel< fem_types >::is_boundary_node ( std::size_t  element_node  )  const [inline]

Returns true if element_node refers to a node index on the reference element which lies at the boundary of the grid in the current element.

Referenced by imaging::SimpleEquationAdaptor< equation_t >::force_vector(), and imaging::SimpleEquationAdaptor< equation_t >::stiffness_matrix().

template<class fem_types>
const ublas::fixed_vector<float_t, fem_types::data_dimension> imaging::FemKernel< fem_types >::boundary_normal_at_node ( std::size_t  element_node  )  const [inline]

Returns the boundary normal in element_node of the current element. The parameter element_node refers to the index of the node as in the reference element.


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

Generated on Tue Feb 10 10:01:31 2009 for imaging2 by  doxygen 1.5.5