Reference documentation for deal.II version 9.1.0-pre
Public Member Functions | Public Attributes | List of all members
MappingQGeneric< dim, spacedim >::InternalData Class Reference

#include <deal.II/fe/mapping_q_generic.h>

Inheritance diagram for MappingQGeneric< dim, spacedim >::InternalData:
[legend]

Public Member Functions

 InternalData (const unsigned int polynomial_degree)
 
void initialize (const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
 
void initialize_face (const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
 
void compute_shape_function_values (const std::vector< Point< dim >> &unit_points)
 
const double & shape (const unsigned int qpoint, const unsigned int shape_nr) const
 
double & shape (const unsigned int qpoint, const unsigned int shape_nr)
 
const Tensor< 1, dim > & derivative (const unsigned int qpoint, const unsigned int shape_nr) const
 
Tensor< 1, dim > & derivative (const unsigned int qpoint, const unsigned int shape_nr)
 
const Tensor< 2, dim > & second_derivative (const unsigned int qpoint, const unsigned int shape_nr) const
 
Tensor< 2, dim > & second_derivative (const unsigned int qpoint, const unsigned int shape_nr)
 
const Tensor< 3, dim > & third_derivative (const unsigned int qpoint, const unsigned int shape_nr) const
 
Tensor< 3, dim > & third_derivative (const unsigned int qpoint, const unsigned int shape_nr)
 
const Tensor< 4, dim > & fourth_derivative (const unsigned int qpoint, const unsigned int shape_nr) const
 
Tensor< 4, dim > & fourth_derivative (const unsigned int qpoint, const unsigned int shape_nr)
 
virtual std::size_t memory_consumption () const override
 
- Public Member Functions inherited from Mapping< dim, spacedim >::InternalDataBase
 InternalDataBase ()
 
 InternalDataBase (const InternalDataBase &)=delete
 
virtual ~InternalDataBase ()=default
 

Public Attributes

std::vector< double > shape_values
 
std::vector< Tensor< 1, dim > > shape_derivatives
 
std::vector< Tensor< 2, dim > > shape_second_derivatives
 
std::vector< Tensor< 3, dim > > shape_third_derivatives
 
std::vector< Tensor< 4, dim > > shape_fourth_derivatives
 
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim-1)> unit_tangentials
 
const unsigned int polynomial_degree
 
const unsigned int n_shape_functions
 
internal::MatrixFreeFunctions::ShapeInfo< VectorizedArray< double > > shape_info
 
AlignedVector< VectorizedArray< double > > scratch
 
AlignedVector< VectorizedArray< double > > values_dofs
 
AlignedVector< VectorizedArray< double > > values_quad
 
AlignedVector< VectorizedArray< double > > gradients_quad
 
AlignedVector< VectorizedArray< double > > hessians_quad
 
bool tensor_product_quadrature
 
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
 
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
 
std::vector< std::vector< Tensor< 1, spacedim > > > aux
 
std::vector< Point< spacedim > > mapping_support_points
 
Triangulation< dim, spacedim >::cell_iterator cell_of_current_support_points
 
std::vector< double > volume_elements
 
- Public Attributes inherited from Mapping< dim, spacedim >::InternalDataBase
UpdateFlags update_each
 

Detailed Description

template<int dim, int spacedim = dim>
class MappingQGeneric< dim, spacedim >::InternalData

Storage for internal data of polynomial mappings. See Mapping::InternalDataBase for an extensive description.

For the current class, the InternalData class stores data that is computed once when the object is created (in get_data()) as well as data the class wants to store from between the call to fill_fe_values(), fill_fe_face_values(), or fill_fe_subface_values() until possible later calls from the finite element to functions such as transform(). The latter class of member variables are marked as 'mutable'.

Definition at line 248 of file mapping_q_generic.h.

Constructor & Destructor Documentation

template<int dim, int spacedim>
MappingQGeneric< dim, spacedim >::InternalData::InternalData ( const unsigned int  polynomial_degree)

Constructor. The argument denotes the polynomial degree of the mapping to which this object will correspond.

Definition at line 650 of file mapping_q_generic.cc.

Member Function Documentation

template<int dim, int spacedim>
void MappingQGeneric< dim, spacedim >::InternalData::initialize ( const UpdateFlags  update_flags,
const Quadrature< dim > &  quadrature,
const unsigned int  n_original_q_points 
)

Initialize the object's member variables related to cell data based on the given arguments.

The function also calls compute_shape_function_values() to actually set the member variables related to the values and derivatives of the mapping shape functions.

Definition at line 682 of file mapping_q_generic.cc.

template<int dim, int spacedim>
void MappingQGeneric< dim, spacedim >::InternalData::initialize_face ( const UpdateFlags  update_flags,
const Quadrature< dim > &  quadrature,
const unsigned int  n_original_q_points 
)

Initialize the object's member variables related to cell and face data based on the given arguments. In order to initialize cell data, this function calls initialize().

Definition at line 795 of file mapping_q_generic.cc.

template<int dim, int spacedim>
void MappingQGeneric< dim, spacedim >::InternalData::compute_shape_function_values ( const std::vector< Point< dim >> &  unit_points)

Compute the values and/or derivatives of the shape functions used for the mapping.

Which values, derivatives, or higher order derivatives are computed is determined by which of the member arrays have nonzero sizes. They are typically set to their appropriate sizes by the initialize() and initialize_face() functions, which indeed call this function internally. However, it is possible (and at times useful) to do the resizing by hand and then call this function directly. An example is in a Newton iteration where we update the location of a quadrature point (e.g., in MappingQ::transform_real_to_uni_cell()) and need to re- compute the mapping and its derivatives at this location, but have already sized all internal arrays correctly.

Definition at line 955 of file mapping_q_generic.cc.

template<int dim, int spacedim = dim>
const double& MappingQGeneric< dim, spacedim >::InternalData::shape ( const unsigned int  qpoint,
const unsigned int  shape_nr 
) const

Shape function at quadrature point. Shape functions are in tensor product order, so vertices must be reordered to obtain transformation.

template<int dim, int spacedim = dim>
double& MappingQGeneric< dim, spacedim >::InternalData::shape ( const unsigned int  qpoint,
const unsigned int  shape_nr 
)

Shape function at quadrature point. See above.

template<int dim, int spacedim = dim>
const Tensor<1, dim>& MappingQGeneric< dim, spacedim >::InternalData::derivative ( const unsigned int  qpoint,
const unsigned int  shape_nr 
) const

Gradient of shape function in quadrature point. See above.

template<int dim, int spacedim = dim>
Tensor<1, dim>& MappingQGeneric< dim, spacedim >::InternalData::derivative ( const unsigned int  qpoint,
const unsigned int  shape_nr 
)

Gradient of shape function in quadrature point. See above.

template<int dim, int spacedim = dim>
const Tensor<2, dim>& MappingQGeneric< dim, spacedim >::InternalData::second_derivative ( const unsigned int  qpoint,
const unsigned int  shape_nr 
) const

Second derivative of shape function in quadrature point. See above.

template<int dim, int spacedim = dim>
Tensor<2, dim>& MappingQGeneric< dim, spacedim >::InternalData::second_derivative ( const unsigned int  qpoint,
const unsigned int  shape_nr 
)

Second derivative of shape function in quadrature point. See above.

template<int dim, int spacedim = dim>
const Tensor<3, dim>& MappingQGeneric< dim, spacedim >::InternalData::third_derivative ( const unsigned int  qpoint,
const unsigned int  shape_nr 
) const

third derivative of shape function in quadrature point. See above.

template<int dim, int spacedim = dim>
Tensor<3, dim>& MappingQGeneric< dim, spacedim >::InternalData::third_derivative ( const unsigned int  qpoint,
const unsigned int  shape_nr 
)

third derivative of shape function in quadrature point. See above.

template<int dim, int spacedim = dim>
const Tensor<4, dim>& MappingQGeneric< dim, spacedim >::InternalData::fourth_derivative ( const unsigned int  qpoint,
const unsigned int  shape_nr 
) const

fourth derivative of shape function in quadrature point. See above.

template<int dim, int spacedim = dim>
Tensor<4, dim>& MappingQGeneric< dim, spacedim >::InternalData::fourth_derivative ( const unsigned int  qpoint,
const unsigned int  shape_nr 
)

fourth derivative of shape function in quadrature point. See above.

template<int dim, int spacedim>
std::size_t MappingQGeneric< dim, spacedim >::InternalData::memory_consumption ( ) const
overridevirtual

Return an estimate (in bytes) or the memory consumption of this object.

Reimplemented from Mapping< dim, spacedim >::InternalDataBase.

Definition at line 662 of file mapping_q_generic.cc.

Member Data Documentation

template<int dim, int spacedim = dim>
std::vector<double> MappingQGeneric< dim, spacedim >::InternalData::shape_values

Values of shape functions. Access by function shape.

Computed once.

Definition at line 374 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
std::vector<Tensor<1, dim> > MappingQGeneric< dim, spacedim >::InternalData::shape_derivatives

Values of shape function derivatives. Access by function derivative.

Computed once.

Definition at line 381 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
std::vector<Tensor<2, dim> > MappingQGeneric< dim, spacedim >::InternalData::shape_second_derivatives

Values of shape function second derivatives. Access by function second_derivative.

Computed once.

Definition at line 389 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
std::vector<Tensor<3, dim> > MappingQGeneric< dim, spacedim >::InternalData::shape_third_derivatives

Values of shape function third derivatives. Access by function second_derivative.

Computed once.

Definition at line 397 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
std::vector<Tensor<4, dim> > MappingQGeneric< dim, spacedim >::InternalData::shape_fourth_derivatives

Values of shape function fourth derivatives. Access by function second_derivative.

Computed once.

Definition at line 405 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
std::array<std::vector<Tensor<1, dim> >, GeometryInfo<dim>::faces_per_cell *(dim - 1)> MappingQGeneric< dim, spacedim >::InternalData::unit_tangentials

Unit tangential vectors. Used for the computation of boundary forms and normal vectors.

This array has (dim-1)*GeometryInfo::faces_per_cell entries. The first GeometryInfo::faces_per_cell contain the vectors in the first tangential direction for each face; the second set of GeometryInfo::faces_per_cell entries contain the vectors in the second tangential direction (only in 3d, since there we have 2 tangential directions per face), etc.

Filled once.

Definition at line 422 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
const unsigned int MappingQGeneric< dim, spacedim >::InternalData::polynomial_degree

The polynomial degree of the mapping. Since the objects here are also used (with minor adjustments) by MappingQ, we need to store this.

Definition at line 428 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
const unsigned int MappingQGeneric< dim, spacedim >::InternalData::n_shape_functions

Number of shape functions. If this is a Q1 mapping, then it is simply the number of vertices per cell. However, since also derived classes use this class (e.g. the Mapping_Q() class), the number of shape functions may also be different.

In general, it is \((p+1)^\text{dim}\), where \(p\) is the polynomial degree of the mapping.

Definition at line 439 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
internal::MatrixFreeFunctions::ShapeInfo<VectorizedArray<double> > MappingQGeneric< dim, spacedim >::InternalData::shape_info

In case the quadrature rule given represents a tensor product we need to store the evaluations of the 1d polynomials at the the 1d quadrature points. That is what this variable is for.

Definition at line 457 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
AlignedVector<VectorizedArray<double> > MappingQGeneric< dim, spacedim >::InternalData::scratch
mutable

In case the quadrature rule given represents a tensor product we need to store temporary data in this object.

Definition at line 463 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
AlignedVector<VectorizedArray<double> > MappingQGeneric< dim, spacedim >::InternalData::values_dofs
mutable

In case the quadrature rule given represents a tensor product the values at the mapped support points are stored in this object.

Definition at line 469 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
AlignedVector<VectorizedArray<double> > MappingQGeneric< dim, spacedim >::InternalData::values_quad
mutable

In case the quadrature rule given represents a tensor product the values at the quadrature points are stored in this object.

Definition at line 475 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
AlignedVector<VectorizedArray<double> > MappingQGeneric< dim, spacedim >::InternalData::gradients_quad
mutable

In case the quadrature rule given represents a tensor product the gradients at the quadrature points are stored in this object.

Definition at line 481 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
AlignedVector<VectorizedArray<double> > MappingQGeneric< dim, spacedim >::InternalData::hessians_quad
mutable

In case the quadrature rule given represents a tensor product the hessians at the quadrature points are stored in this object.

Definition at line 487 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
bool MappingQGeneric< dim, spacedim >::InternalData::tensor_product_quadrature

Indicates whether the given Quadrature object is a tensor product.

Definition at line 492 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
std::vector<DerivativeForm<1, dim, spacedim> > MappingQGeneric< dim, spacedim >::InternalData::covariant
mutable

Tensors of covariant transformation at each of the quadrature points. The matrix stored is the Jacobian * G^{-1}, where G = Jacobian^{t} * Jacobian, is the first fundamental form of the map; if dim=spacedim then it reduces to the transpose of the inverse of the Jacobian matrix, which itself is stored in the contravariant field of this structure.

Computed on each cell.

Definition at line 503 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
std::vector<DerivativeForm<1, dim, spacedim> > MappingQGeneric< dim, spacedim >::InternalData::contravariant
mutable

Tensors of contravariant transformation at each of the quadrature points. The contravariant matrix is the Jacobian of the transformation, i.e. \(J_{ij}=dx_i/d\hat x_j\).

Computed on each cell.

Definition at line 512 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
std::vector<std::vector<Tensor<1, spacedim> > > MappingQGeneric< dim, spacedim >::InternalData::aux
mutable

Auxiliary vectors for internal use.

Definition at line 517 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
std::vector<Point<spacedim> > MappingQGeneric< dim, spacedim >::InternalData::mapping_support_points
mutable

Stores the support points of the mapping shape functions on the cell_of_current_support_points.

Definition at line 523 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
Triangulation<dim, spacedim>::cell_iterator MappingQGeneric< dim, spacedim >::InternalData::cell_of_current_support_points
mutable

Stores the cell of which the mapping_support_points are stored.

Definition at line 529 of file mapping_q_generic.h.

template<int dim, int spacedim = dim>
std::vector<double> MappingQGeneric< dim, spacedim >::InternalData::volume_elements
mutable

The determinant of the Jacobian in each quadrature point. Filled if update_volume_elements.

Definition at line 535 of file mapping_q_generic.h.


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