
Public Member Functions | |
| virtual | ~Mapping () |
| virtual Point< dim > | transform_unit_to_real_cell (const typename Triangulation< dim >::cell_iterator &cell, const Point< dim > &p) const =0 |
| virtual Point< dim > | transform_real_to_unit_cell (const typename Triangulation< dim >::cell_iterator &cell, const Point< dim > &p) const =0 |
| virtual void | transform_covariant (const VectorSlice< const std::vector< Tensor< 1, dim > > > input, const unsigned int offset, VectorSlice< std::vector< Tensor< 1, dim > > > output, const InternalDataBase &internal) const =0 |
| virtual void | transform_covariant (const VectorSlice< const std::vector< Tensor< 2, dim > > > input, const unsigned int offset, VectorSlice< std::vector< Tensor< 2, dim > > > output, const InternalDataBase &internal) const =0 |
| virtual void | transform_contravariant (const VectorSlice< const std::vector< Tensor< 1, dim > > > input, const unsigned int offset, VectorSlice< std::vector< Tensor< 1, dim > > > output, const typename Mapping< dim >::InternalDataBase &internal) const =0 |
| virtual void | transform_contravariant (const VectorSlice< const std::vector< Tensor< 2, dim > > > intput, const unsigned int offset, const VectorSlice< std::vector< Tensor< 2, dim > > > output, const typename Mapping< dim >::InternalDataBase &internal) const =0 |
| const Point< dim > & | support_point_value (const unsigned int index, const typename Mapping< dim >::InternalDataBase &internal) const |
| const Tensor< 2, dim > & | support_point_gradient (const unsigned int index, const typename Mapping< dim >::InternalDataBase &internal) const |
| const Tensor< 2, dim > & | support_point_inverse_gradient (const unsigned int index, const typename Mapping< dim >::InternalDataBase &internal) const |
| virtual Mapping< dim > * | clone () const =0 |
Private Member Functions | |
| virtual UpdateFlags | update_once (const UpdateFlags) const =0 |
| virtual UpdateFlags | update_each (const UpdateFlags) const =0 |
| virtual InternalDataBase * | get_data (const UpdateFlags, const Quadrature< dim > &quadrature) const =0 |
| virtual InternalDataBase * | get_face_data (const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const =0 |
| virtual InternalDataBase * | get_subface_data (const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const =0 |
| virtual void | fill_fe_values (const typename Triangulation< dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, InternalDataBase &internal, std::vector< Point< dim > > &quadrature_points, std::vector< double > &JxW_values, std::vector< Tensor< 2, dim > > &jacobians, std::vector< Tensor< 3, dim > > &jacobian_grads, std::vector< Tensor< 2, dim > > &inverse_jacobians) const =0 |
| virtual void | fill_fe_face_values (const typename Triangulation< dim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, InternalDataBase &internal, std::vector< Point< dim > > &quadrature_points, std::vector< double > &JxW_values, std::vector< Tensor< 1, dim > > &boundary_form, std::vector< Point< dim > > &normal_vectors, std::vector< double > &cell_JxW_values) const =0 |
| virtual void | fill_fe_subface_values (const typename Triangulation< dim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim-1 > &quadrature, InternalDataBase &internal, std::vector< Point< dim > > &quadrature_points, std::vector< double > &JxW_values, std::vector< Tensor< 1, dim > > &boundary_form, std::vector< Point< dim > > &normal_vectors, std::vector< double > &cell_JxW_values) const =0 |
Friends | |
| class | FEValuesBase< dim > |
| class | FEValues< dim > |
| class | FEFaceValues< dim > |
| class | FESubfaceValues< dim > |
Classes | |
| class | ExcInvalidData |
| class | InternalDataBase |
The interface for filling the tables of FEValues is provided. Everything else has to happen in derived classes.
The following paragraph applies to the implementation of FEValues. Usage of the class is as follows: first, call the functions update_once and update_each with the update flags you need. This includes the flags needed by the FiniteElement. Then call get_*_data and with the or'd results. This will initialize and return some internal data structures. On the first cell, call fill_fe_*_values with the result of update_once. Finally, on each cell, use fill_fe_*_values with the result of update_each to compute values for a special cell.
A hint to implementators: no function except the two functions update_once and update_each may add any flags.
| virtual Point<dim> Mapping< dim >::transform_unit_to_real_cell | ( | const typename Triangulation< dim >::cell_iterator & | cell, | |
| const Point< dim > & | p | |||
| ) | const [pure virtual] |
Transforms the point p on the unit cell to the point p_real on the real cell cell and returns p_real.
Implemented in MappingCartesian< dim >, MappingQ< dim >, and MappingQ1< dim >.
| virtual Point<dim> Mapping< dim >::transform_real_to_unit_cell | ( | const typename Triangulation< dim >::cell_iterator & | cell, | |
| const Point< dim > & | p | |||
| ) | const [pure virtual] |
Transforms the point p on the real cell to the point p_unit on the unit cell cell and returns p_unit.
Implemented in MappingCartesian< dim >, MappingQ< dim >, and MappingQ1< dim >.
| virtual void Mapping< dim >::transform_covariant | ( | const VectorSlice< const std::vector< Tensor< 1, dim > > > | input, | |
| const unsigned int | offset, | |||
| VectorSlice< std::vector< Tensor< 1, dim > > > | output, | |||
| const InternalDataBase & | internal | |||
| ) | const [pure virtual] |
Transform a field of covariant vectors. The covariant transformation multiplies a vector from the right with the inverse of the Jacobian matrix of the transformation from unit to real space cell. Alternatively, this is equivalent to a multiplication from the left with the transpose if the inverse matrix.
The list of arguments is as follows: we transform as many elements in the input field, starting from offset, as there are elements in output. The input array may hold more elements than are needed (some finite element classes use this for a denser storage of data), but it needs to have at least output.size() elements starting with element offset.
| virtual void Mapping< dim >::transform_covariant | ( | const VectorSlice< const std::vector< Tensor< 2, dim > > > | input, | |
| const unsigned int | offset, | |||
| VectorSlice< std::vector< Tensor< 2, dim > > > | output, | |||
| const InternalDataBase & | internal | |||
| ) | const [pure virtual] |
Transform a set of matrices covariantly, i.e. multiply each function in the input range by the Jacobian matrices at the different quadrature points from the left.
The same applies as to the function above regarding input and output arguments.
| virtual void Mapping< dim >::transform_contravariant | ( | const VectorSlice< const std::vector< Tensor< 1, dim > > > | input, | |
| const unsigned int | offset, | |||
| VectorSlice< std::vector< Tensor< 1, dim > > > | output, | |||
| const typename Mapping< dim >::InternalDataBase & | internal | |||
| ) | const [pure virtual] |
Transform a field of contravariant vectors. The contravariant transformation multiplies a vector from the left with the Jacobian matrix of the transformation from unit to real space cell.
The list of arguments is as follows: we transform as many elements in the input field, starting from offset, as there are elements in output. The input array may hold more elements than are needed (some finite element classes use this for a denser storage of data), but it needs to have at least output.size() elements starting with element offset.
Implemented in MappingCartesian< dim >, MappingQ< dim >, and MappingQ1< dim >.
| virtual void Mapping< dim >::transform_contravariant | ( | const VectorSlice< const std::vector< Tensor< 2, dim > > > | intput, | |
| const unsigned int | offset, | |||
| const VectorSlice< std::vector< Tensor< 2, dim > > > | output, | |||
| const typename Mapping< dim >::InternalDataBase & | internal | |||
| ) | const [pure virtual] |
Transform a set of matrices contravariantly, i.e. multiply each function in the input range by the inverse Jacobian matrices at the different quadrature points from the right. Note that here it is no more equivalent to multiplication with the transpose of the inverse matrices from the left, unless the matrices to be multiplied with are symmetric.
The same applies as to the function above regarding input and output arguments.
Implemented in MappingCartesian< dim >, MappingQ< dim >, and MappingQ1< dim >.
| const Point< dim > & Mapping< dim >::support_point_value | ( | const unsigned int | index, | |
| const typename Mapping< dim >::InternalDataBase & | internal | |||
| ) | const [inline] |
The transformed (generalized) support point.
References AssertIndexRange.
| const Tensor< 2, dim > & Mapping< dim >::support_point_gradient | ( | const unsigned int | index, | |
| const typename Mapping< dim >::InternalDataBase & | internal | |||
| ) | const [inline] |
The Jacobian matrix of the transformation in the (generalized) support point.
References AssertIndexRange.
| const Tensor< 2, dim > & Mapping< dim >::support_point_inverse_gradient | ( | const unsigned int | index, | |
| const typename Mapping< dim >::InternalDataBase & | internal | |||
| ) | const [inline] |
The inverse Jacobian matrix of the transformation in the (generalized) support point.
References AssertIndexRange.
Return a pointer to a copy of the present object. The caller of this copy then assumes ownership of it.
Since one can't create objects of class Mapping, this function of course has to be implemented by derived classes.
This function is mainly used by the hp::MappingCollection class.
Implemented in MappingC1< dim >, MappingCartesian< dim >, MappingQ< dim >, MappingQ1< dim >, MappingQ1Eulerian< dim, EulerVectorType >, and MappingQEulerian< dim, EulerVectorType >.
| virtual UpdateFlags Mapping< dim >::update_once | ( | const | UpdateFlags | ) | const [private, pure 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.
Implemented in MappingCartesian< dim >, and MappingQ1< dim >.
| virtual UpdateFlags Mapping< dim >::update_each | ( | const | UpdateFlags | ) | const [private, pure 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.
Implemented in MappingCartesian< dim >, and MappingQ1< dim >.
| virtual InternalDataBase* Mapping< dim >::get_data | ( | const | UpdateFlags, | |
| const Quadrature< dim > & | quadrature | |||
| ) | const [private, pure virtual] |
Prepare internal data structures and fill in values independent of the cell.
Implemented in MappingCartesian< dim >, MappingQ< dim >, and MappingQ1< dim >.
| virtual InternalDataBase* Mapping< dim >::get_face_data | ( | const UpdateFlags | flags, | |
| const Quadrature< dim-1 > & | quadrature | |||
| ) | const [private, pure virtual] |
Prepare internal data structure for transformation of faces and fill in values independent of the cell.
Implemented in MappingCartesian< dim >, MappingQ< dim >, and MappingQ1< dim >.
| virtual InternalDataBase* Mapping< dim >::get_subface_data | ( | const UpdateFlags | flags, | |
| const Quadrature< dim-1 > & | quadrature | |||
| ) | const [private, pure virtual] |
Prepare internal data structure for transformation of children of faces and fill in values independent of the cell.
Implemented in MappingCartesian< dim >, MappingQ< dim >, and MappingQ1< dim >.
| virtual void Mapping< dim >::fill_fe_values | ( | const typename Triangulation< dim >::cell_iterator & | cell, | |
| const Quadrature< dim > & | quadrature, | |||
| InternalDataBase & | internal, | |||
| std::vector< Point< dim > > & | quadrature_points, | |||
| std::vector< double > & | JxW_values, | |||
| std::vector< Tensor< 2, dim > > & | jacobians, | |||
| std::vector< Tensor< 3, dim > > & | jacobian_grads, | |||
| std::vector< Tensor< 2, dim > > & | inverse_jacobians | |||
| ) | const [private, pure virtual] |
Fill the transformation fields of FEValues. Given a grid cell and the quadrature points on the unit cell, it computes all values specified by flags. The arrays to be filled have to have the correct size.
Values are split into two groups: first, quadrature_points and JxW_values are filled with the quadrature rule transformed to the cell in physical space.
The second group contains the matrices needed to transform vector-valued functions, namely jacobians, the derivatives jacobian_grads, and the inverse operations in inverse_jacobians.
| virtual void Mapping< dim >::fill_fe_face_values | ( | const typename Triangulation< dim >::cell_iterator & | cell, | |
| const unsigned int | face_no, | |||
| const Quadrature< dim-1 > & | quadrature, | |||
| InternalDataBase & | internal, | |||
| std::vector< Point< dim > > & | quadrature_points, | |||
| std::vector< double > & | JxW_values, | |||
| std::vector< Tensor< 1, dim > > & | boundary_form, | |||
| std::vector< Point< dim > > & | normal_vectors, | |||
| std::vector< double > & | cell_JxW_values | |||
| ) | const [private, pure virtual] |
Performs the same as fill_fe_values on a face. Additionally, boundary_form and normal_vectors can be computed on surfaces. The boundary form is the vector product of the image of coordinate vectors on the surface of the unit cell. It is a vector normal to the surface, pointing outwards and having the length of the surface element. Therefore, it is more economic to use the boundary form instead of the product of the unit normal and the transformed quadrature weight.
| virtual void Mapping< dim >::fill_fe_subface_values | ( | const typename Triangulation< dim >::cell_iterator & | cell, | |
| const unsigned int | face_no, | |||
| const unsigned int | sub_no, | |||
| const Quadrature< dim-1 > & | quadrature, | |||
| InternalDataBase & | internal, | |||
| std::vector< Point< dim > > & | quadrature_points, | |||
| std::vector< double > & | JxW_values, | |||
| std::vector< Tensor< 1, dim > > & | boundary_form, | |||
| std::vector< Point< dim > > & | normal_vectors, | |||
| std::vector< double > & | cell_JxW_values | |||
| ) | const [private, pure virtual] |
See above.
friend class FEValuesBase< dim > [friend] |
Give class FEValues access to the private get_...data and fill_fe_...values functions.
friend class FEValues< dim > [friend] |
friend class FEFaceValues< dim > [friend] |
friend class FESubfaceValues< dim > [friend] |
documentation generated on Fri Nov 21 07:13:15 2008 by
doxygen
1.5.6