Classes | |
class | InternalData |
Public Member Functions | |
virtual Mapping< dim, spacedim > ::InternalDataBase * | get_data (const UpdateFlags, const Quadrature< dim > &quadrature) const |
virtual Mapping< dim, spacedim > ::InternalDataBase * | get_face_data (const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const |
virtual Mapping< dim, spacedim > ::InternalDataBase * | get_subface_data (const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const |
virtual void | fill_fe_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Quadrature< dim > &quadrature, typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, std::vector< Point< spacedim > > &quadrature_points, std::vector< double > &JxW_values, std::vector< Tensor< 2, spacedim > > &jacobians, std::vector< Tensor< 3, spacedim > > &jacobian_grads, std::vector< Tensor< 2, spacedim > > &inverse_jacobians, std::vector< Point< spacedim > > &, enum CellSimilarity::Similarity &cell_similarity) const |
virtual void | fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, std::vector< Point< dim > > &quadrature_points, std::vector< double > &JxW_values, std::vector< Tensor< 1, dim > > &boundary_form, std::vector< Point< spacedim > > &normal_vectors) const |
virtual void | fill_fe_subface_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim-1 > &quadrature, typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, std::vector< Point< dim > > &quadrature_points, std::vector< double > &JxW_values, std::vector< Tensor< 1, dim > > &boundary_form, std::vector< Point< spacedim > > &normal_vectors) const |
virtual void | transform (const VectorSlice< const std::vector< Tensor< 1, dim > > > input, VectorSlice< std::vector< Tensor< 1, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const MappingType type) const |
virtual void | transform (const VectorSlice< const std::vector< Tensor< 2, dim > > > input, VectorSlice< std::vector< Tensor< 2, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const MappingType type) const |
virtual Point< spacedim > | transform_unit_to_real_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const |
virtual Point< dim > | transform_real_to_unit_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const |
virtual Mapping< dim, spacedim > * | clone () const |
Protected Member Functions | |
void | compute_fill (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const enum CellSimilarity::Similarity cell_similarity, InternalData &data, std::vector< Point< dim > > &quadrature_points, std::vector< Point< dim > > &normal_vectors) const |
Private Member Functions | |
virtual UpdateFlags | update_once (const UpdateFlags) const |
virtual UpdateFlags | update_each (const UpdateFlags) const |
Static Private Attributes | |
static const unsigned int | invalid_face_number = numbers::invalid_unsigned_int |
This class maps the unit cell to a grid cell with surfaces parallel to the coordinate lines/planes. The mapping is therefore a scaling along the coordinate directions. It is specifically developed for cartesian meshes. Apply this mapping to a general mesh to get strange results.
For more information about the spacedim
template parameter check the documentation of FiniteElement or the one of Triangulation.
virtual Mapping<dim, spacedim>::InternalDataBase* MappingCartesian< dim, spacedim >::get_data | ( | const | UpdateFlags, | |
const Quadrature< dim > & | quadrature | |||
) | const [virtual] |
Prepare internal data structures and fill in values independent of the cell.
Implements Mapping< dim, spacedim >.
virtual Mapping<dim, spacedim>::InternalDataBase* MappingCartesian< dim, spacedim >::get_face_data | ( | const UpdateFlags | flags, | |
const Quadrature< dim-1 > & | quadrature | |||
) | const [virtual] |
Prepare internal data structure for transformation of faces and fill in values independent of the cell.
Implements Mapping< dim, spacedim >.
virtual Mapping<dim, spacedim>::InternalDataBase* MappingCartesian< dim, spacedim >::get_subface_data | ( | const UpdateFlags | flags, | |
const Quadrature< dim-1 > & | quadrature | |||
) | const [virtual] |
Prepare internal data structure for transformation of children of faces and fill in values independent of the cell.
Implements Mapping< dim, spacedim >.
virtual void MappingCartesian< dim, spacedim >::fill_fe_values | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, | |
const Quadrature< dim > & | quadrature, | |||
typename Mapping< dim, spacedim >::InternalDataBase & | mapping_data, | |||
std::vector< Point< spacedim > > & | quadrature_points, | |||
std::vector< double > & | JxW_values, | |||
std::vector< Tensor< 2, spacedim > > & | jacobians, | |||
std::vector< Tensor< 3, spacedim > > & | jacobian_grads, | |||
std::vector< Tensor< 2, spacedim > > & | inverse_jacobians, | |||
std::vector< Point< spacedim > > & | , | |||
enum CellSimilarity::Similarity & | cell_similarity | |||
) | const [virtual] |
virtual void MappingCartesian< dim, spacedim >::fill_fe_face_values | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, | |
const unsigned int | face_no, | |||
const Quadrature< dim-1 > & | quadrature, | |||
typename Mapping< dim, spacedim >::InternalDataBase & | mapping_data, | |||
std::vector< Point< dim > > & | quadrature_points, | |||
std::vector< double > & | JxW_values, | |||
std::vector< Tensor< 1, dim > > & | boundary_form, | |||
std::vector< Point< spacedim > > & | normal_vectors | |||
) | const [virtual] |
virtual void MappingCartesian< dim, spacedim >::fill_fe_subface_values | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, | |
const unsigned int | face_no, | |||
const unsigned int | sub_no, | |||
const Quadrature< dim-1 > & | quadrature, | |||
typename Mapping< dim, spacedim >::InternalDataBase & | mapping_data, | |||
std::vector< Point< dim > > & | quadrature_points, | |||
std::vector< double > & | JxW_values, | |||
std::vector< Tensor< 1, dim > > & | boundary_form, | |||
std::vector< Point< spacedim > > & | normal_vectors | |||
) | const [virtual] |
virtual void MappingCartesian< dim, spacedim >::transform | ( | const VectorSlice< const std::vector< Tensor< 1, dim > > > | input, | |
VectorSlice< std::vector< Tensor< 1, spacedim > > > | output, | |||
const typename Mapping< dim, spacedim >::InternalDataBase & | internal, | |||
const MappingType | type | |||
) | const [virtual] |
virtual void MappingCartesian< dim, spacedim >::transform | ( | const VectorSlice< const std::vector< Tensor< 2, dim > > > | input, | |
VectorSlice< std::vector< Tensor< 2, spacedim > > > | output, | |||
const typename Mapping< dim, spacedim >::InternalDataBase & | internal, | |||
const MappingType | type | |||
) | const [virtual] |
virtual Point<spacedim> MappingCartesian< dim, spacedim >::transform_unit_to_real_cell | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, | |
const Point< dim > & | p | |||
) | const [virtual] |
Transforms the point p
on the unit cell to the point p_real
on the real cell cell
and returns p_real
.
Implements Mapping< dim, spacedim >.
virtual Point<dim> MappingCartesian< dim, spacedim >::transform_real_to_unit_cell | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, | |
const Point< spacedim > & | p | |||
) | const [virtual] |
Transforms the point p
on the real cell to the point p_unit
on the unit cell cell
and returns p_unit
.
Uses Newton iteration and the transform_unit_to_real_cell
function.
Implements Mapping< dim, spacedim >.
virtual Mapping<dim, spacedim>* MappingCartesian< dim, spacedim >::clone | ( | ) | const [virtual] |
Return a pointer to a copy of the present object. The caller of this copy then assumes ownership of it.
Implements Mapping< dim, spacedim >.
void MappingCartesian< dim, spacedim >::compute_fill | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, | |
const unsigned int | face_no, | |||
const unsigned int | sub_no, | |||
const enum CellSimilarity::Similarity | cell_similarity, | |||
InternalData & | data, | |||
std::vector< Point< dim > > & | quadrature_points, | |||
std::vector< Point< dim > > & | normal_vectors | |||
) | const [protected] |
Do the computation for the fill_*
functions.
virtual UpdateFlags MappingCartesian< dim, spacedim >::update_once | ( | const | UpdateFlags | ) | const [private, virtual] |
Indicate fields to be updated in the constructor of FEValues. Especially, fields not asked for by FEValues, but computed for efficiency reasons will be notified here.
See The interplay of UpdateFlags, Mapping and FiniteElement in FEValues.
Implements Mapping< dim, spacedim >.
virtual UpdateFlags MappingCartesian< dim, spacedim >::update_each | ( | const | UpdateFlags | ) | const [private, virtual] |
The same as update_once(), but for the flags to be updated for each grid cell.
See The interplay of UpdateFlags, Mapping and FiniteElement in FEValues.
Implements Mapping< dim, spacedim >.
const unsigned int MappingCartesian< dim, spacedim >::invalid_face_number = numbers::invalid_unsigned_int [static, private] |
Value to indicate that a given face or subface number is invalid.