Reference documentation for deal.II version 8.5.1
Public Member Functions | Private Attributes | List of all members
internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number > Class Template Reference

#include <deal.II/matrix_free/mapping_data_on_the_fly.h>

Public Member Functions

 MappingDataOnTheFly (const Mapping< dim > &mapping, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags)
 
 MappingDataOnTheFly (const Quadrature< 1 > &quadrature, const UpdateFlags update_flags)
 
void reinit (typename ::Triangulation< dim >::cell_iterator cell)
 
bool is_initialized () const
 
typename ::Triangulation< dim >::cell_iterator get_cell () const
 
const ::FEValues< dim > & get_fe_values () const
 
const AlignedVector< Tensor< 2, dim, VectorizedArray< Number > > > & get_inverse_jacobians () const
 
const AlignedVector< VectorizedArray< Number > > & get_JxW_values () const
 
const AlignedVector< Point< dim, VectorizedArray< Number > > > & get_quadrature_points () const
 
const AlignedVector< Tensor< 1, dim, VectorizedArray< Number > > > & get_normal_vectors () const
 
const Quadrature< 1 > & get_quadrature () const
 

Private Attributes

typename ::Triangulation< dim >::cell_iterator present_cell
 
FE_Nothing< dim > fe_dummy
 
::FEValues< dim > fe_values
 
const Quadrature< 1 > quadrature_1d
 
AlignedVector< Tensor< 2, dim, VectorizedArray< Number > > > inverse_jacobians
 
AlignedVector< VectorizedArray< Number > > jxw_values
 
AlignedVector< Point< dim, VectorizedArray< Number > > > quadrature_points
 
AlignedVector< Tensor< 1, dim, VectorizedArray< Number > > > normal_vectors
 

Detailed Description

template<int dim, typename Number = double>
class internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >

This class provides evaluated mapping information using standard deal.II information in a form that FEEvaluation and friends can use for vectorized access. Since no vectorization over cells is available with the DoFHandler/Triangulation cell iterators, the interface to FEEvaluation's vectorization model is to use VectorizedArray::n_array_element copies of the same element. This interface is thus primarily useful for evaluating several operators on the same cell, e.g., when assembling cell matrices.

As opposed to the Mapping classes in deal.II, this class does not actually provide a boundary description that can be used to evaluate the geometry, but it rather provides the evaluated geometry from a given deal.II mapping (as passed to the constructor of this class) in a form accessible to FEEvaluation.

Author
Martin Kronbichler, 2014

Definition at line 58 of file mapping_data_on_the_fly.h.

Constructor & Destructor Documentation

◆ MappingDataOnTheFly() [1/2]

template<int dim, typename Number >
internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::MappingDataOnTheFly ( const Mapping< dim > &  mapping,
const Quadrature< 1 > &  quadrature,
const UpdateFlags  update_flags 
)
inline

Constructor, similar to FEValues. Since this class only evaluates the geometry, no finite element has to be specified and the simplest element, FE_Nothing, is used internally for the underlying FEValues object.

Definition at line 190 of file mapping_data_on_the_fly.h.

◆ MappingDataOnTheFly() [2/2]

template<int dim, typename Number >
internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::MappingDataOnTheFly ( const Quadrature< 1 > &  quadrature,
const UpdateFlags  update_flags 
)
inline

Constructor. This constructor is equivalent to the other one except that it makes the object use a \(Q_1\) mapping (i.e., an object of type MappingQGeneric(1)) implicitly.

Definition at line 210 of file mapping_data_on_the_fly.h.

Member Function Documentation

◆ reinit()

template<int dim, typename Number >
void internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::reinit ( typename ::Triangulation< dim >::cell_iterator  cell)
inline

Initialize with the given cell iterator.

Definition at line 230 of file mapping_data_on_the_fly.h.

◆ is_initialized()

template<int dim, typename Number >
bool internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::is_initialized ( ) const
inline

Return whether reinit() has been called at least once, i.e., a cell has been set.

Definition at line 258 of file mapping_data_on_the_fly.h.

◆ get_cell()

template<int dim, typename Number >
typename::Triangulation< dim >::cell_iterator internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::get_cell ( ) const
inline

Return a triangulation iterator to the current cell.

Definition at line 268 of file mapping_data_on_the_fly.h.

◆ get_fe_values()

template<int dim, typename Number >
const ::FEValues< dim > & internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::get_fe_values ( ) const
inline

Return a reference to the underlying FEValues object that evaluates certain quantities (only mapping-related ones like Jacobians or mapped quadrature points are accessible, as no finite element data is actually used).

Definition at line 278 of file mapping_data_on_the_fly.h.

◆ get_inverse_jacobians()

template<int dim, typename Number >
const AlignedVector< Tensor< 2, dim, VectorizedArray< Number > > > & internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::get_inverse_jacobians ( ) const
inline

Return a vector of inverse transpose Jacobians. For compatibility with FEEvaluation, it returns tensors of vectorized arrays, even though all components are equal.

Definition at line 288 of file mapping_data_on_the_fly.h.

◆ get_JxW_values()

template<int dim, typename Number >
const AlignedVector< VectorizedArray< Number > > & internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::get_JxW_values ( ) const
inline

Return a vector of quadrature weights times the Jacobian determinant (JxW). For compatibility with FEEvaluation, it returns tensors of vectorized arrays, even though all components are equal.

Definition at line 318 of file mapping_data_on_the_fly.h.

◆ get_quadrature_points()

template<int dim, typename Number >
const AlignedVector< Point< dim, VectorizedArray< Number > > > & internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::get_quadrature_points ( ) const
inline

Return a vector of quadrature points in real space on the given cell. For compatibility with FEEvaluation, it returns tensors of vectorized arrays, even though all components are equal.

Definition at line 308 of file mapping_data_on_the_fly.h.

◆ get_normal_vectors()

template<int dim, typename Number >
const AlignedVector< Tensor< 1, dim, VectorizedArray< Number > > > & internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::get_normal_vectors ( ) const
inline

Return a vector of normal vectors in real space on the given cell. For compatibility with FEEvaluation, it returns tensors of vectorized arrays, even though all components are equal.

Definition at line 298 of file mapping_data_on_the_fly.h.

◆ get_quadrature()

template<int dim, typename Number >
const Quadrature< 1 > & internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::get_quadrature ( ) const
inline

Return a reference to 1D quadrature underlying this object.

Definition at line 328 of file mapping_data_on_the_fly.h.

Member Data Documentation

◆ present_cell

template<int dim, typename Number = double>
typename ::Triangulation<dim>::cell_iterator internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::present_cell
private

A cell iterator in case we generate the data on the fly to be able to check if we need to re-generate the information stored in this class.

Definition at line 146 of file mapping_data_on_the_fly.h.

◆ fe_dummy

template<int dim, typename Number = double>
FE_Nothing<dim> internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::fe_dummy
private

Dummy finite element object necessary for initializing the FEValues object.

Definition at line 152 of file mapping_data_on_the_fly.h.

◆ fe_values

template<int dim, typename Number = double>
::FEValues<dim> internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::fe_values
private

An underlying FEValues object that performs the (scalar) evaluation.

Definition at line 157 of file mapping_data_on_the_fly.h.

◆ quadrature_1d

template<int dim, typename Number = double>
const Quadrature<1> internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::quadrature_1d
private

Get 1D quadrature formula to be used for reinitializing shape info.

Definition at line 162 of file mapping_data_on_the_fly.h.

◆ inverse_jacobians

template<int dim, typename Number = double>
AlignedVector<Tensor<2,dim,VectorizedArray<Number> > > internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::inverse_jacobians
private

Inverse Jacobians, stored in vectorized array form.

Definition at line 167 of file mapping_data_on_the_fly.h.

◆ jxw_values

template<int dim, typename Number = double>
AlignedVector<VectorizedArray<Number> > internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::jxw_values
private

Stored Jacobian determinants and quadrature weights

Definition at line 172 of file mapping_data_on_the_fly.h.

◆ quadrature_points

template<int dim, typename Number = double>
AlignedVector<Point<dim,VectorizedArray<Number> > > internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::quadrature_points
private

Stored quadrature points

Definition at line 177 of file mapping_data_on_the_fly.h.

◆ normal_vectors

template<int dim, typename Number = double>
AlignedVector<Tensor<1,dim,VectorizedArray<Number> > > internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number >::normal_vectors
private

Stored normal vectors (for face integration)

Definition at line 182 of file mapping_data_on_the_fly.h.


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