![]() |
Reference documentation for deal.II version 8.1.0
|
#include <mapping_q1.h>
Classes | |
| class | InternalData |
Public Types | |
| typedef QProjector< dim > ::DataSetDescriptor | DataSetDescriptor |
Public Member Functions | |
| MappingQ1 () | |
| 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 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< DerivativeForm< 1, dim, spacedim > > > input, VectorSlice< std::vector< Tensor< 2, 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 Mapping< dim, spacedim > * | clone () 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, typename std::vector< Point< spacedim > > &quadrature_points, std::vector< double > &JxW_values, std::vector< DerivativeForm< 1, dim, spacedim > > &jacobians, std::vector< DerivativeForm< 2, dim, spacedim > > &jacobian_grads, std::vector< DerivativeForm< 1, spacedim, dim > > &inverse_jacobians, std::vector< Point< spacedim > > &cell_normal_vectors, 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, typename std::vector< Point< spacedim > > &quadrature_points, std::vector< double > &JxW_values, typename std::vector< Tensor< 1, spacedim > > &boundary_form, typename 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, typename std::vector< Point< spacedim > > &quadrature_points, std::vector< double > &JxW_values, typename std::vector< Tensor< 1, spacedim > > &boundary_form, typename std::vector< Point< spacedim > > &normal_vectors) const |
| void | compute_shapes (const std::vector< Point< dim > > &unit_points, InternalData &data) const |
| void | compute_data (const UpdateFlags flags, const Quadrature< dim > &quadrature, const unsigned int n_orig_q_points, InternalData &data) const |
| void | compute_face_data (const UpdateFlags flags, const Quadrature< dim > &quadrature, const unsigned int n_orig_q_points, InternalData &data) const |
| void | compute_fill (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int npts, const DataSetDescriptor data_set, const CellSimilarity::Similarity cell_similarity, InternalData &data, std::vector< Point< spacedim > > &quadrature_points) const |
| void | compute_fill_face (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const unsigned int npts, const DataSetDescriptor data_set, const std::vector< double > &weights, InternalData &mapping_data, std::vector< Point< spacedim > > &quadrature_points, std::vector< double > &JxW_values, std::vector< Tensor< 1, spacedim > > &boundary_form, std::vector< Point< spacedim > > &normal_vectors) const |
| virtual void | compute_shapes_virtual (const std::vector< Point< dim > > &unit_points, InternalData &data) const |
| Point< spacedim > | transform_unit_to_real_cell_internal (const InternalData &mdata) const |
| Point< dim > | transform_real_to_unit_cell_internal (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit, InternalData &mdata) const |
| virtual bool | preserves_vertex_locations () const |
| template<> | |
| Point< 2 > | transform_real_to_unit_cell_internal (const Triangulation< 2, 3 >::cell_iterator &cell, const Point< 3 > &p, const Point< 2 > &initial_p_unit, InternalData &mdata) const |
| template<> | |
| Point< 1 > | transform_real_to_unit_cell_internal (const Triangulation< 1, 2 >::cell_iterator &cell, const Point< 2 > &p, const Point< 1 > &initial_p_unit, InternalData &mdata) const |
| template<> | |
| Point< 1 > | transform_real_to_unit_cell_internal (const Triangulation< 1, 3 >::cell_iterator &cell, const Point< 3 > &p, const Point< 1 > &initial_p_unit, InternalData &mdata) const |
Public Member Functions inherited from Mapping< dim, spacedim > | |
| virtual | ~Mapping () |
| virtual void | transform (const VectorSlice< const std::vector< Tensor< 1, dim > > > input, VectorSlice< std::vector< Tensor< 1, spacedim > > > output, const InternalDataBase &internal, const MappingType type) const =0 |
| virtual void | transform (const VectorSlice< const std::vector< DerivativeForm< 1, dim, spacedim > > > input, VectorSlice< std::vector< Tensor< 2, spacedim > > > output, const InternalDataBase &internal, const MappingType type) const =0 |
| virtual void | transform (const VectorSlice< const std::vector< Tensor< 2, dim > > > input, VectorSlice< std::vector< Tensor< 2, spacedim > > > output, const InternalDataBase &internal, const MappingType type) const =0 |
| void | transform_covariant (const VectorSlice< const std::vector< Tensor< 1, dim > > > input, const unsigned int offset, VectorSlice< std::vector< Tensor< 1, spacedim > > > output, const InternalDataBase &internal) const DEAL_II_DEPRECATED |
| void | transform_covariant (const VectorSlice< const std::vector< DerivativeForm< 1, dim, spacedim > > > input, const unsigned int offset, VectorSlice< std::vector< Tensor< 2, spacedim > > > output, const InternalDataBase &internal) const DEAL_II_DEPRECATED |
| void | transform_contravariant (const VectorSlice< const std::vector< Tensor< 1, dim > > > input, const unsigned int offset, VectorSlice< std::vector< Tensor< 1, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal) const DEAL_II_DEPRECATED |
| void | transform_contravariant (const VectorSlice< const std::vector< DerivativeForm< 1, dim, spacedim > > > input, const unsigned int offset, const VectorSlice< std::vector< Tensor< 2, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal) const DEAL_II_DEPRECATED |
| const Point< spacedim > & | support_point_value (const unsigned int index, const typename Mapping< dim, spacedim >::InternalDataBase &internal) const |
| const Tensor< 2, spacedim > & | support_point_gradient (const unsigned int index, const typename Mapping< dim, spacedim >::InternalDataBase &internal) const |
| const Tensor< 2, spacedim > & | support_point_inverse_gradient (const unsigned int index, const typename Mapping< dim, spacedim >::InternalDataBase &internal) const |
| DeclException0 (ExcInvalidData) | |
| DeclException0 (ExcTransformationFailed) | |
| DeclException3 (ExcDistortedMappedCell, Point< spacedim >, double, int,<< "The image of the mapping applied to cell with center ["<< arg1<< "] is distorted. The cell geometry or the "<< "mapping are invalid, giving a non-positive volume "<< "fraction of "<< arg2<< " in quadrature point "<< arg3<< ".") | |
Public Member Functions inherited from Subscriptor | |
| Subscriptor () | |
| Subscriptor (const Subscriptor &) | |
| virtual | ~Subscriptor () |
| Subscriptor & | operator= (const Subscriptor &) |
| void | subscribe (const char *identifier=0) const |
| void | unsubscribe (const char *identifier=0) const |
| unsigned int | n_subscriptions () const |
| void | list_subscribers () const |
| DeclException3 (ExcInUse, int, char *, std::string &,<< "Object of class "<< arg2<< " is still used by "<< arg1<< " other objects.\n"<< "(Additional information: "<< arg3<< ")\n"<< "Note the entry in the Frequently Asked Questions of "<< "deal.II (linked to from http://www.dealii.org/) for "<< "more information on what this error means.") | |
| DeclException2 (ExcNoSubscriber, char *, char *,<< "No subscriber with identifier \""<< arg2<< "\" did subscribe to this object of class "<< arg1) | |
| template<class Archive > | |
| void | serialize (Archive &ar, const unsigned int version) |
Protected Member Functions | |
| template<int rank> | |
| void | transform_fields (const VectorSlice< const std::vector< Tensor< rank, dim > > > input, VectorSlice< std::vector< Tensor< rank, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const MappingType type) const |
| template<int rank> | |
| void | transform_gradients (const VectorSlice< const std::vector< Tensor< rank, dim > > > input, VectorSlice< std::vector< Tensor< rank, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const MappingType type) const |
| template<int rank> | |
| void | transform_differential_forms (const VectorSlice< const std::vector< DerivativeForm< rank, dim, spacedim > > > input, VectorSlice< std::vector< DerivativeForm< rank, spacedim, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const MappingType type) const |
| template<int dim_> | |
| Point< dim_ > | transform_real_to_unit_cell_internal_codim1 (const typename Triangulation< dim_, dim_+1 >::cell_iterator &cell, const Point< dim_+1 > &p, const Point< dim_ > &initial_p_unit, InternalData &mdata) const |
| Point< dim > | transform_real_to_unit_cell_initial_guess (const std::vector< Point< spacedim > > &vertex, const Point< spacedim > &p) const |
Private Member Functions | |
| virtual UpdateFlags | update_once (const UpdateFlags flags) const |
| virtual UpdateFlags | update_each (const UpdateFlags flags) const |
| 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 | compute_mapping_support_points (const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const |
Static Private Attributes | |
| static const unsigned int | n_shape_functions = GeometryInfo<dim>::vertices_per_cell |
Mapping of general quadrilateral/hexahedra by d-linear shape functions.
This function maps the unit cell to a general grid cell with straight lines in
dimensions (remark that in 3D the surfaces may be curved, even if the edges are not). This is the well-known mapping for polyhedral domains.
Shape function for this mapping are the same as for the finite element FE_Q of order 1. Therefore, coupling these two yields an isoparametric element.
For more information about the spacedim template parameter check the documentation of FiniteElement or the one of Triangulation.
Definition at line 57 of file mapping_q1.h.
| typedef QProjector<dim>::DataSetDescriptor MappingQ1< dim, spacedim >::DataSetDescriptor |
Declare a convenience typedef for the class that describes offsets into quadrature formulas projected onto faces and subfaces.
Definition at line 391 of file mapping_q1.h.
Default constructor.
|
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 >.
Reimplemented in MappingQ< dim, spacedim >.
|
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.
In the codimension one case, this function returns the normal projection of the real point p on the curve or surface identified by the cell.
p. If this is the case then this function throws an exception of type Mapping::ExcTransformationFailed . Whether the given point p lies outside the cell can therefore be determined by checking whether the return reference coordinates lie inside of outside the reference cell (e.g., using GeometryInfo::is_inside_unit_cell) or whether the exception mentioned above has been thrown. Implements Mapping< dim, spacedim >.
Reimplemented in MappingQ< dim, spacedim >.
|
protected |
This function and the next allow to generate the transform require by the virtual transform() in mapping, but unfortunately in C++ one cannot declare a virtual template function.
|
protected |
see doc in transform_fields
|
protected |
see doc in transform_fields
|
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 >.
Reimplemented in MappingQ< dim, spacedim >, MappingQ1Eulerian< dim, VECTOR, spacedim >, MappingQEulerian< dim, VECTOR, spacedim >, and MappingC1< dim, spacedim >.
|
virtual |
Implementation of the interface in Mapping.
Reimplemented in MappingQ< dim, spacedim >, MappingQ1Eulerian< dim, VECTOR, spacedim >, and MappingQEulerian< dim, VECTOR, spacedim >.
|
virtual |
Implementation of the interface in Mapping.
Reimplemented in MappingQ< dim, spacedim >.
|
virtual |
Implementation of the interface in Mapping.
Reimplemented in MappingQ< dim, spacedim >.
| void MappingQ1< dim, spacedim >::compute_shapes | ( | const std::vector< Point< dim > > & | unit_points, |
| InternalData & | data | ||
| ) | const |
Compute shape values and/or derivatives.
Calls either the compute_shapes_virtual of this class or that of the derived class, depending on whether data.is_mapping_q1_data equals true or false.
| void MappingQ1< dim, spacedim >::compute_data | ( | const UpdateFlags | flags, |
| const Quadrature< dim > & | quadrature, | ||
| const unsigned int | n_orig_q_points, | ||
| InternalData & | data | ||
| ) | const |
Do the computations for the get_data functions. Here, the data vectors of InternalData are reinitialized to proper size and shape values are computed.
| void MappingQ1< dim, spacedim >::compute_face_data | ( | const UpdateFlags | flags, |
| const Quadrature< dim > & | quadrature, | ||
| const unsigned int | n_orig_q_points, | ||
| InternalData & | data | ||
| ) | const |
Do the computations for the get_face_data functions. Here, the data vectors of InternalData are reinitialized to proper size and shape values and derivatives are computed. Furthermore unit_tangential vectors of the face are computed.
| void MappingQ1< dim, spacedim >::compute_fill | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, |
| const unsigned int | npts, | ||
| const DataSetDescriptor | data_set, | ||
| const CellSimilarity::Similarity | cell_similarity, | ||
| InternalData & | data, | ||
| std::vector< Point< spacedim > > & | quadrature_points | ||
| ) | const |
Do the computation for the fill_* functions.
| void MappingQ1< dim, spacedim >::compute_fill_face | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, |
| const unsigned int | face_no, | ||
| const unsigned int | subface_no, | ||
| const unsigned int | npts, | ||
| const DataSetDescriptor | data_set, | ||
| const std::vector< double > & | weights, | ||
| InternalData & | mapping_data, | ||
| std::vector< Point< spacedim > > & | quadrature_points, | ||
| std::vector< double > & | JxW_values, | ||
| std::vector< Tensor< 1, spacedim > > & | boundary_form, | ||
| std::vector< Point< spacedim > > & | normal_vectors | ||
| ) | const |
Do the computation for the fill_* functions.
|
virtual |
Compute shape values and/or derivatives.
| Point<spacedim> MappingQ1< dim, spacedim >::transform_unit_to_real_cell_internal | ( | const InternalData & | mdata | ) | const |
Transforms a point p on the unit cell to the point p_real on the real cell cell and returns p_real.
This function is called by transform_unit_to_real_cell and multiple times (through the Newton iteration) by transform_real_to_unit_cell_internal.
Takes a reference to an InternalData that must already include the shape values at point p and the mapping support points of the cell.
This InternalData argument avoids multiple computations of the shape values at point p and especially multiple computations of the mapping support points.
| Point<dim> MappingQ1< dim, spacedim >::transform_real_to_unit_cell_internal | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, |
| const Point< spacedim > & | p, | ||
| const Point< dim > & | initial_p_unit, | ||
| InternalData & | mdata | ||
| ) | const |
Transforms the point p on the real cell to the corresponding point on the unit cell cell by a Newton iteration.
Takes a reference to an InternalData that is assumed to be previously created by the get_data function with UpdateFlags including update_transformation_values and update_transformation_gradients and a one point Quadrature that includes the given initial guess for the transformation initial_p_unit. Hence this function assumes that mdata already includes the transformation shape values and gradients computed at initial_p_unit.
mdata will be changed by this function.
|
virtual |
Always returns true because MappingQ1 preserves vertex locations.
Implements Mapping< dim, spacedim >.
Reimplemented in MappingQ1Eulerian< dim, VECTOR, spacedim >, and MappingQEulerian< dim, VECTOR, spacedim >.
|
protected |
Compute an initial guess to pass to the Newton method in transform_real_to_unit_cell.
For the initial guess we proceed in the following way:
Some details about how we compute the least square plane. We look for a spacedim x (dim + 1) matrix X such that
X * M = Y
where M is a (dim+1) x n_vertices matrix and Y a spacedim x n_vertices.
And: The i-th column of M is unit_vertex[i] and the last row all 1's. The i-th column of Y is real_vertex[i].
If we split X=[A|b], the least square approx is A x_hat+b
Classically X = Y * (M^t (M M^t)^{-1})
Let K = M^t * (M M^t)^{-1} = [KA Kb] this can be precomputed, and that is exactely what we do.
Finally A = Y*KA and b = Y*Kb.
|
privatevirtual |
Implementation of the interface in Mapping.
Description of effects:
update_quadrature_points is required, the output will contain update_transformation_values. This computes the values of the transformation basis polynomials at the unit cell quadrature points. update_covariant_transformation, update_contravariant_transformation, update_JxW_values, update_boundary_forms, update_normal_vectors is required, the output will contain update_transformation_gradients to compute derivatives of the transformation basis polynomials. Implements Mapping< dim, spacedim >.
|
privatevirtual |
Implementation of the interface in Mapping.
Description of effects if flags contains:
update_quadrature_points is copied to the output to compute the quadrature points on the real cell. update_JxW_values is copied and requires update_boundary_forms on faces. The latter, because the surface element is just the norm of the boundary form. update_normal_vectors is copied and requires update_boundary_forms. The latter, because the normal vector is the normalized boundary form. update_covariant_transformation is copied and requires update_contravariant_transformation, since it is computed as the inverse of the latter. update_JxW_values is copied and requires update_contravariant_transformation, since it is computed as one over determinant of the latter. update_boundary_forms is copied and requires update_contravariant_transformation, since the boundary form is computed as the contravariant image of the normal vector to the unit cell. Implements Mapping< dim, spacedim >.
|
privatevirtual |
Prepare internal data structures and fill in values independent of the cell.
Implements Mapping< dim, spacedim >.
Reimplemented in MappingQ< dim, spacedim >.
|
privatevirtual |
Prepare internal data structure for transformation of faces and fill in values independent of the cell.
Implements Mapping< dim, spacedim >.
Reimplemented in MappingQ< dim, spacedim >.
|
privatevirtual |
Prepare internal data structure for transformation of children of faces and fill in values independent of the cell.
Implements Mapping< dim, spacedim >.
Reimplemented in MappingQ< dim, spacedim >.
|
privatevirtual |
Computes the support points of the mapping. For MappingQ1 these are the vertices. However, other classes may override this function. In particular, the MappingQ1Eulerian class does exactly this by not computing the support points from the geometry of the current cell but instead evaluating an externally given displacement field in addition to the geometry of the cell.
Reimplemented in MappingQ< dim, spacedim >, MappingQEulerian< dim, VECTOR, spacedim >, and MappingQ1Eulerian< dim, VECTOR, spacedim >.
|
staticprivate |
Number of shape functions. Is simply the number of vertices per cell for the Q1 mapping.
Definition at line 760 of file mapping_q1.h.
1.8.6