#include <Grid.hpp>
Public Types | |
typedef ublas::fixed_vector < float_t, fem_types::data_dimension > | vertex_t |
The data type of the element coordinates (determined by the template parameter fem_types). | |
typedef fem_types::transform_t | transform_t |
The element transformation class (determined by the template parameter fem_types). | |
typedef fem_types::shape_function_t | shape_function_t |
The shape function class (determined by the template parameter fem_types). | |
typedef fem_types::integrator_t | integrator_t |
The integrator class (determined by the template parameter fem_types). | |
typedef fem_types::boundary_integrator_t | boundary_integrator_t |
The boundary integrator class (determined by the template parameter fem_types). | |
Public Member Functions | |
Grid () | |
size_t | n_elements () const |
size_t | n_vertices () const |
size_t | n_nodes () const |
size_t | n_boundary_nodes () const |
size_t | n_boundary_elements () const |
void | set_vertex (size_t global_vertex_index, const vertex_t &vertex) |
const vertex_t & | vertex (size_t global_vertex_index) const |
void | set_vertex (size_t element_index, size_t vertex_index, const vertex_t &vertex) |
const vertex_t & | vertex (size_t element_index, size_t vertex_index) const |
size_t | global_vertex_index (size_t element_index, size_t vertex_index) const |
void | set_global_vertex_index (size_t element_index, size_t vertex_index, size_t global_vertex_index) |
size_t | global_node_index (size_t element_index, size_t node_index) const |
void | set_global_node_index (size_t element_index, size_t node_index, size_t global_node_index) |
void | set_boundary_element (size_t boundary_element_index, size_t parent_element_index, size_t parent_element_face) |
size_t | parent_element (size_t boundary_element_index) const |
size_t | parent_element_face (size_t boundary_element_index) const |
bool | is_boundary_node (size_t global_node_index) const |
bool | is_boundary_node (size_t element_index, size_t node_index) const |
const vertex_t | boundary_normal (size_t global_node_index) const |
const vertex_t | boundary_normal (size_t element_index, size_t node_index) const |
void | set_regular (bool is_regular) |
bool | is_regular () const |
void | set_dimensions (size_t n_vertices, size_t n_elements, size_t n_boundary_elements, size_t n_nodes) |
void | element_transform (size_t element_index, transform_t &transform) const |
void | interpolate_value (size_t integrator_node, const ublas::vector< float_t > &data, const FemKernel< fem_types > &kernel, float_t &value) const |
template<size_t N> | |
void | interpolate_value (size_t integrator_node, const ublas::vector< float_t > &data, const FemKernel< fem_types > &kernel, ublas::fixed_vector< float_t, N > &value) const |
void | interpolate_gradient (size_t integrator_node, const ublas::vector< float_t > &data, const FemKernel< fem_types > &kernel, ublas::fixed_vector< float_t, data_dimension > &gradient) const |
template<size_t N> | |
void | interpolate_gradient (size_t integrator_node, const ublas::vector< float_t > &data, const FemKernel< fem_types > &kernel, ublas::fixed_matrix< float_t, data_dimension, N > &gradient) const |
void | interpolate_boundary_value (size_t boundary_integrator_node, const ublas::mapped_vector< float_t > &data, const FemKernel< fem_types > &kernel, float_t &value) const |
template<size_t N> | |
void | interpolate_boundary_value (size_t boundary_integrator_node, const ublas::mapped_vector< float_t > &data, const FemKernel< fem_types > &kernel, ublas::fixed_vector< float_t, N > &value) const |
void | interpolate_boundary_derivative (size_t boundary_integrator_node, const ublas::mapped_vector< float_t > &data, const FemKernel< fem_types > &kernel, float_t &value) const |
template<size_t N> | |
void | interpolate_boundary_derivative (size_t boundary_integrator_node, const ublas::mapped_vector< float_t > &data, const FemKernel< fem_types > &kernel, ublas::fixed_vector< float_t, N > &value) const |
void | global_coordinates (size_t integrator_node, const FemKernel< fem_types > &kernel, ublas::fixed_vector< float_t, data_dimension > &coordinates) const |
Static Public Attributes | |
static const size_t | data_dimension = fem_types::data_dimension |
The number of coordinates of the element vertices (determined by the template parameter fem_types). | |
static const size_t | n_element_vertices = fem_types::transform_t::n_element_vertices |
The number of vertices per element (determined by the template parameter fem_types). | |
static const size_t | n_element_nodes = fem_types::shape_function_t::n_element_nodes |
The number of nodes per element (determined by the template parameter fem_types). | |
static const size_t | n_element_faces = fem_types::transform_t::n_element_faces |
The number of nodes per element (determined by the template parameter fem_types). | |
Classes | |
class | BoundaryElement |
Grid stores all the data to access the elements of a FE grid and obtain their geometry. In addition it keeps track of the boundary of the grid and elements there.
In the following we give a more detailed description of how Grid stores this information. The basic building block of a grid is an element. We assume an element to have vertices, faces and nodes. Each of these has a unique index (running from zero to Grid::n_element_vertices, from zero to Grid::n_element_faces and from zero to Grid::n_element_nodes) on the reference element which carries over to the transformed elements in the grid.
An element in the grid is identified by its element index which runs from zero to n_elements(). Every element has Grid::n_element_vertices vertices which are identified by their vertex indices. Each of the vertex indices is mapped to a global vertex index which refers to its position in a list of vertex coordinates.In other words, a combination of an element index + vertex index can be mapped to the global vertex index which again identifies the position of the element vertex coordinates in the vertex coordinate list. The vertex coordinates determine the position and the geometry of the element. In a grid consisting of triangular elements an element will always have three vertices. Elements in a quadrilateral or tetrahedral grid are defined by 4 vertices respectively. This also means that the number of vertices per element and the dimension of the vertex coordinates (which are determined by Grid::n_element_vertices and Grid::data_dimension) always have to correspond to the element transformation (Grid::transform_t).
The element nodes, in contrast, are defined by the shape function on the element. Again each node of a element is determined by its node index and is further mapped to a global node index. The coordinates of the nodes do not have to stored because the are completely determined by the coordinates of the vertices of their element. Instead, the global node indices refer to the position of the node in the stiffness matrix and the force vector of the FE problem. For linear shape functions on triangle elements or bilinear shape functions on quadrilateral elements the number of nodes per element is the same as the number of vertices and their actual positions coincide. For higher order shape functions there maybe more nodes than vertices. Clearly, the number of nodes per element (determined by Grid::n_element_nodes) has to match the type of shape functions (Grid::shape_function_t).
Finally, Grid keeps track of the boundary of the FE domain. The term boundary element refers to a face of an element (the parent element) which lies a the boundary of the grid. It is identified by its element index which runs from zero to n_boundary_elements(). For each boundary element the index of its parent element and the index of the corresponding face on the boundary element is stored. A boundary vertex is an element vertex which happens to on a boundary element. The same holds for a boundary node. For each boundary node the outer unit normal at the grid boundary is stored. This information has to be provided during the construction of the grid.
imaging::Grid< fem_types >::Grid | ( | ) | [inline] |
Default constructor.
size_t imaging::Grid< fem_types >::n_elements | ( | ) | const [inline] |
Returns the number of elements of the grid (excluding boundary elements).
Referenced by imaging::Assembler::assemble(), imaging::Assembler::assemble_force_vector(), imaging::Assembler::assemble_stiffness_matrix(), imaging::average(), imaging::integrate(), imaging::maximum(), imaging::minimum(), and imaging::triangulate_shape().
size_t imaging::Grid< fem_types >::n_vertices | ( | ) | const [inline] |
Returns the number of element vertices of the grid. In case there is exactly shape function one per element vertex this equals the number of nodes.
Referenced by imaging::triangulate_shape().
size_t imaging::Grid< fem_types >::n_nodes | ( | ) | const [inline] |
Returns the number of element nodes of this grid. The number of element nodes equals the number of degrees of freedom of FE discretization of a scalar PDE on this grid. For a system of equations this number has to be multiplied by the number of equations to obtain the total number of degrees of freedom.
Referenced by imaging::Assembler::assemble(), imaging::Assembler::assemble_force_vector(), and imaging::Assembler::assemble_stiffness_matrix().
size_t imaging::Grid< fem_types >::n_boundary_nodes | ( | ) | const [inline] |
Returns the number of nodes which are at the boundary of the grid. The equations for these nodes will be influenced by boundary conditions.
size_t imaging::Grid< fem_types >::n_boundary_elements | ( | ) | const [inline] |
Returns the number of boundary elements, i.e. those boundary elements of parent elements which are at the actual boundary of the grid.
Referenced by imaging::Assembler::assemble(), imaging::Assembler::assemble_force_vector(), and imaging::Assembler::assemble_stiffness_matrix().
void imaging::Grid< fem_types >::set_vertex | ( | size_t | global_vertex_index, | |
const vertex_t & | vertex | |||
) | [inline] |
Sets the vertex with global index global_vertex_index to vertex.
Referenced by imaging::ellipse_grid(), imaging::triangulate_shape(), and imaging::uniform_grid().
const vertex_t& imaging::Grid< fem_types >::vertex | ( | size_t | global_vertex_index | ) | const [inline] |
Returns a const reference to the vertex with global index global_vertex_index.
void imaging::Grid< fem_types >::set_vertex | ( | size_t | element_index, | |
size_t | vertex_index, | |||
const vertex_t & | vertex | |||
) | [inline] |
Sets the vertex with index vertex_index on the element element_index to vertex.
References imaging::Grid< fem_types >::global_vertex_index().
const vertex_t& imaging::Grid< fem_types >::vertex | ( | size_t | element_index, | |
size_t | vertex_index | |||
) | const [inline] |
Returns a const reference to the vertex with index vertex_index on the element element_index.
References imaging::Grid< fem_types >::global_vertex_index().
size_t imaging::Grid< fem_types >::global_vertex_index | ( | size_t | element_index, | |
size_t | vertex_index | |||
) | const [inline] |
Returns the global vertex index of the node with index vertex_index on the element element_index. This index corresponds to the index of the vertex coordinates in the vertex coordinate list stored by the Grid object.
Referenced by imaging::Grid< fem_types >::element_transform(), imaging::Grid< fem_types >::set_vertex(), and imaging::Grid< fem_types >::vertex().
void imaging::Grid< fem_types >::set_global_vertex_index | ( | size_t | element_index, | |
size_t | vertex_index, | |||
size_t | global_vertex_index | |||
) | [inline] |
Sets the global vertex index of the vertex with index vertex_index on the element element_index. This index corresponds to the index of the vertex coordinates in the vertex coordinate list stored by the Grid object.
Referenced by imaging::ellipse_grid(), imaging::triangulate_shape(), and imaging::uniform_grid().
size_t imaging::Grid< fem_types >::global_node_index | ( | size_t | element_index, | |
size_t | node_index | |||
) | const [inline] |
Returns the global node index of the node with index node_index on the element element_index. This index corresponds to the position of the node in the stiffness matrix and the force vector of the associated FE problem.
Referenced by imaging::Assembler::assemble(), imaging::Assembler::assemble_force_vector(), imaging::Assembler::assemble_stiffness_matrix(), imaging::Grid< fem_types >::boundary_normal(), imaging::Grid< fem_types >::interpolate_boundary_derivative(), imaging::Grid< fem_types >::interpolate_boundary_value(), imaging::Grid< fem_types >::interpolate_gradient(), imaging::Grid< fem_types >::interpolate_value(), imaging::Grid< fem_types >::is_boundary_node(), and imaging::Grid< fem_types >::set_boundary_element().
void imaging::Grid< fem_types >::set_global_node_index | ( | size_t | element_index, | |
size_t | node_index, | |||
size_t | global_node_index | |||
) | [inline] |
Sets the global node index of the node with index node_index on the element element_index. This index corresponds to the position of the node in the stiffness matrix and the force vector of the associated FE problem.
Referenced by imaging::ellipse_grid(), imaging::triangulate_shape(), and imaging::uniform_grid().
void imaging::Grid< fem_types >::set_boundary_element | ( | size_t | boundary_element_index, | |
size_t | parent_element_index, | |||
size_t | parent_element_face | |||
) | [inline] |
Sets the boundary element boundary_element_index. Note that the nodes on the boundary element are not automatically set to be boundary nodes; the user has to do this manually during grid construction.
References imaging::Grid< fem_types >::element_transform(), and imaging::Grid< fem_types >::global_node_index().
Referenced by imaging::ellipse_grid(), imaging::triangulate_shape(), and imaging::uniform_grid().
size_t imaging::Grid< fem_types >::parent_element | ( | size_t | boundary_element_index | ) | const [inline] |
Returns the index of the parent element of the boundary element boundary_element_index.
Referenced by imaging::Assembler::assemble(), imaging::Assembler::assemble_force_vector(), imaging::Assembler::assemble_stiffness_matrix(), imaging::Grid< fem_types >::interpolate_boundary_derivative(), and imaging::Grid< fem_types >::interpolate_boundary_value().
size_t imaging::Grid< fem_types >::parent_element_face | ( | size_t | boundary_element_index | ) | const [inline] |
Returns the index of the face of the parent element of the boundary element boundary_element_index.
bool imaging::Grid< fem_types >::is_boundary_node | ( | size_t | global_node_index | ) | const [inline] |
Returns true if the node global_node_index lies at the boundary of the grid.
Referenced by imaging::Grid< fem_types >::is_boundary_node().
bool imaging::Grid< fem_types >::is_boundary_node | ( | size_t | element_index, | |
size_t | node_index | |||
) | const [inline] |
Returns true if the node node_index on the element element_index lies at the boundary of the grid.
References imaging::Grid< fem_types >::global_node_index(), and imaging::Grid< fem_types >::is_boundary_node().
const vertex_t imaging::Grid< fem_types >::boundary_normal | ( | size_t | global_node_index | ) | const [inline] |
Returns the unit boundary normal of the node global_node_index. If the node happens not to be a boundary node an Exception is thrown.
Referenced by imaging::Grid< fem_types >::boundary_normal().
const vertex_t imaging::Grid< fem_types >::boundary_normal | ( | size_t | element_index, | |
size_t | node_index | |||
) | const [inline] |
Returns the unit boundary normal of the node node_index on the element element_index. If the node happens not to be a boundary node an Exception is thrown.
References imaging::Grid< fem_types >::boundary_normal(), and imaging::Grid< fem_types >::global_node_index().
void imaging::Grid< fem_types >::set_regular | ( | bool | is_regular | ) | [inline] |
Returns true if the grid is regular, i.e. the geometry of each of its elements is the identical modulo rigid transformation. Marking a grid as regular speeds up the assembly of the stiffness matrix and the force vector. Prominent examples of regular grids are pixel and voxel discretizations.
bool imaging::Grid< fem_types >::is_regular | ( | ) | const [inline] |
Returns true if the grid is regular, i.e. the geometry of each of its elements is the identical modulo rigid transformation. Marking a grid as regular speeds up the assembly of the stiffness matrix and the force vector. Prominent examples of regular grids are pixel and voxel discretizations.
Referenced by imaging::Assembler::assemble(), imaging::Assembler::assemble_force_vector(), imaging::Assembler::assemble_stiffness_matrix(), imaging::average(), imaging::integrate(), imaging::maximum(), and imaging::minimum().
void imaging::Grid< fem_types >::set_dimensions | ( | size_t | n_vertices, | |
size_t | n_elements, | |||
size_t | n_boundary_elements, | |||
size_t | n_nodes | |||
) | [inline] |
Sets the dimensions of the grid. The parameters refer to the total numbers of vertices, elements, boundary elements and nodes, respectively. Note that vertices or nodes which belong to more than one element are only counted once.
Referenced by imaging::ellipse_grid(), imaging::triangulate_shape(), and imaging::uniform_grid().
void imaging::Grid< fem_types >::element_transform | ( | size_t | element_index, | |
transform_t & | transform | |||
) | const [inline] |
Initializes the element transformation transform to the element element_index.
References imaging::Grid< fem_types >::global_vertex_index(), and imaging::Grid< fem_types >::n_element_vertices.
Referenced by imaging::Grid< fem_types >::global_coordinates(), and imaging::Grid< fem_types >::set_boundary_element().
void imaging::Grid< fem_types >::interpolate_value | ( | size_t | integrator_node, | |
const ublas::vector< float_t > & | data, | |||
const FemKernel< fem_types > & | kernel, | |||
float_t & | value | |||
) | const [inline] |
Interpolates data in integrator_node on the current element of kernel. The result is written to value.
[in] | data | A vector of size n_nodes(). Its values can be interpreted as the evaluation of a scalar function on the nodes of the grid. interpolate_value() interpolates these values and evaluates the interpolation. |
[in] | kernel | A FemKernel object. It must be initialized to a valid element of the grid. |
References imaging::FemKernel< fem_types >::current_element(), imaging::Grid< fem_types >::global_node_index(), imaging::Grid< fem_types >::n_element_nodes, and imaging::FemKernel< fem_types >::shape_value().
void imaging::Grid< fem_types >::interpolate_value | ( | size_t | integrator_node, | |
const ublas::vector< float_t > & | data, | |||
const FemKernel< fem_types > & | kernel, | |||
ublas::fixed_vector< float_t, N > & | value | |||
) | const [inline] |
Interpretes data as a vector of N-dimensional vectors and interpolates it in integrator_node on the current element of kernel. The result is written to value.
[in] | data | A vector of size (N * n_nodes()). Its values can be interpreted as the evaluation of an N-dimensional function on the nodes of the grid. interpolate_value() interpolates these values and evaluates the interpolation. |
[in] | kernel | A FemKernel object. It must be initialized to a valid element of the grid. |
References imaging::FemKernel< fem_types >::current_element(), imaging::Grid< fem_types >::global_node_index(), imaging::Grid< fem_types >::n_element_nodes, and imaging::FemKernel< fem_types >::shape_value().
void imaging::Grid< fem_types >::interpolate_gradient | ( | size_t | integrator_node, | |
const ublas::vector< float_t > & | data, | |||
const FemKernel< fem_types > & | kernel, | |||
ublas::fixed_vector< float_t, data_dimension > & | gradient | |||
) | const [inline] |
Interpolates the gradient of data in integrator_node on the current element of kernel. The result is written to value.
[in] | data | A vector of size n_nodes(). Its values can be interpreted as the evaluation of a scalar function on the nodes of the grid. interpolate_value() interpolates these values and evaluates the gradient of the interpolation. |
[in] | kernel | A FemKernel object. It must be initialized to a valid element of the grid. |
References imaging::FemKernel< fem_types >::current_element(), imaging::Grid< fem_types >::global_node_index(), imaging::Grid< fem_types >::n_element_nodes, and imaging::FemKernel< fem_types >::shape_gradient().
void imaging::Grid< fem_types >::interpolate_gradient | ( | size_t | integrator_node, | |
const ublas::vector< float_t > & | data, | |||
const FemKernel< fem_types > & | kernel, | |||
ublas::fixed_matrix< float_t, data_dimension, N > & | gradient | |||
) | const [inline] |
Interpretes data as a vector of N-dimensional vectors and interpolates its gradient in integrator_node on the current element of kernel. The result is written to value.
[in] | data | A vector of size (N * n_nodes()). Its values can be interpreted as the evaluation of an N-dimensional function on the nodes of the grid. interpolate_value() interpolates these values and evaluates the gradient (a matrix) of the interpolation. |
[in] | kernel | A FemKernel object. It must be initialized to a valid element of the grid. |
References imaging::FemKernel< fem_types >::current_element(), imaging::Grid< fem_types >::global_node_index(), imaging::Grid< fem_types >::n_element_nodes, and imaging::FemKernel< fem_types >::shape_gradient().
void imaging::Grid< fem_types >::interpolate_boundary_value | ( | size_t | boundary_integrator_node, | |
const ublas::mapped_vector< float_t > & | data, | |||
const FemKernel< fem_types > & | kernel, | |||
float_t & | value | |||
) | const [inline] |
Interpolates data in boundary_integrator_node on the current boundary element of kernel. The result is written to value.
[in] | data | A sparse vector of size n_nodes(). It must contain values at positions which correspond to boundary nodes. Its values can be interpreted as the evaluation of a scalar function on the boundary nodes of the grid. interpolate_value() interpolates these values and evaluates the interpolation. |
[in] | kernel | A FemKernel object. It must be initialized to a valid boundary element of the grid. |
References imaging::FemKernel< fem_types >::current_boundary_element(), imaging::Grid< fem_types >::global_node_index(), imaging::Grid< fem_types >::n_element_nodes, imaging::Grid< fem_types >::parent_element(), and imaging::FemKernel< fem_types >::shape_boundary_value().
void imaging::Grid< fem_types >::interpolate_boundary_value | ( | size_t | boundary_integrator_node, | |
const ublas::mapped_vector< float_t > & | data, | |||
const FemKernel< fem_types > & | kernel, | |||
ublas::fixed_vector< float_t, N > & | value | |||
) | const [inline] |
Interpretes data as a vector of N-dimensional vectors and interpolates it in boundary_integrator_node on the current boundary element of kernel. The result is written to value.
[in] | data | A sparse vector of size (N * n_nodes()). It must contain values at positions which correspond to boundary nodes. Its values can be interpreted as the evaluation of an N-dimensional function on the boundary nodes of the grid. interpolate_value() interpolates these values and evaluates the interpolation. |
[in] | kernel | A FemKernel object. It must be initialized to a valid boundary element of the grid. |
References imaging::FemKernel< fem_types >::current_boundary_element(), imaging::Grid< fem_types >::global_node_index(), imaging::Grid< fem_types >::n_element_nodes, imaging::Grid< fem_types >::parent_element(), and imaging::FemKernel< fem_types >::shape_boundary_value().
void imaging::Grid< fem_types >::interpolate_boundary_derivative | ( | size_t | boundary_integrator_node, | |
const ublas::mapped_vector< float_t > & | data, | |||
const FemKernel< fem_types > & | kernel, | |||
float_t & | value | |||
) | const [inline] |
Interpolates the tangential derivative of data in boundary_integrator_node on the current boundary element of kernel. The result is written to value.
[in] | data | A sparse vector of size n_nodes(). It must contain values at positions which correspond to boundary nodes. Its values can be interpreted as the evaluation of a scalar function on the boundary nodes of the grid. interpolate_value() interpolates these values and evaluates the tangential derivative of the interpolation. Note that it is not possible to compute the normal derivate if data is known only on the boundary. |
[in] | kernel | A FemKernel object. It must be initialized to a valid boundary element of the grid. |
References imaging::FemKernel< fem_types >::current_boundary_element(), imaging::Grid< fem_types >::global_node_index(), imaging::Grid< fem_types >::n_element_nodes, imaging::Grid< fem_types >::parent_element(), and imaging::FemKernel< fem_types >::shape_boundary_derivative().
void imaging::Grid< fem_types >::interpolate_boundary_derivative | ( | size_t | boundary_integrator_node, | |
const ublas::mapped_vector< float_t > & | data, | |||
const FemKernel< fem_types > & | kernel, | |||
ublas::fixed_vector< float_t, N > & | value | |||
) | const [inline] |
Interpretes data as a vector of N-dimensional vectors and interpolates it in boundary_integrator_node on the current boundary element of kernel. The result is written to value.
[in] | data | A sparse vector of size (N * n_nodes()). It must contain values at positions which correspond to boundary nodes. Its values can be interpreted as the evaluation of a N-dimensional function on the boundary nodes of the grid. interpolate_value() interpolates these values and evaluates the tangential derivative of the interpolation. Note that it is not possible to compute the normal derivate if data is known only on the boundary. |
[in] | kernel | A FemKernel object. It must be initialized to a valid boundary element of the grid. |
References imaging::FemKernel< fem_types >::current_boundary_element(), imaging::Grid< fem_types >::global_node_index(), imaging::Grid< fem_types >::n_element_nodes, imaging::Grid< fem_types >::parent_element(), and imaging::FemKernel< fem_types >::shape_boundary_derivative().
void imaging::Grid< fem_types >::global_coordinates | ( | size_t | integrator_node, | |
const FemKernel< fem_types > & | kernel, | |||
ublas::fixed_vector< float_t, data_dimension > & | coordinates | |||
) | const [inline] |
Computes the global coordinates of the integrator_node on the current element of kernel. The result is written to coordinates.
References imaging::FemKernel< fem_types >::current_element(), and imaging::Grid< fem_types >::element_transform().