Reference documentation for deal.II version 8.5.1
Public Member Functions | Protected Member Functions | List of all members
FEEvaluationAccess< dim, dim, Number > Class Template Reference

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

Inheritance diagram for FEEvaluationAccess< dim, dim, Number >:
[legend]

Public Member Functions

gradient_type get_gradient (const unsigned int q_point) const
 
VectorizedArray< Number > get_divergence (const unsigned int q_point) const
 
SymmetricTensor< 2, dim, VectorizedArray< Number > > get_symmetric_gradient (const unsigned int q_point) const
 
Tensor< 1,(dim==2?1:dim), VectorizedArray< Number > > get_curl (const unsigned int q_point) const
 
Tensor< 3, dim, VectorizedArray< Number > > get_hessian (const unsigned int q_point) const
 
gradient_type get_hessian_diagonal (const unsigned int q_point) const
 
void submit_gradient (const gradient_type grad_in, const unsigned int q_point)
 
void submit_gradient (const Tensor< 1, dim, Tensor< 1, dim, VectorizedArray< Number > > > grad_in, const unsigned int q_point)
 
void submit_divergence (const VectorizedArray< Number > div_in, const unsigned int q_point)
 
void submit_symmetric_gradient (const SymmetricTensor< 2, dim, VectorizedArray< Number > > grad_in, const unsigned int q_point)
 
void submit_curl (const Tensor< 1, dim==2?1:dim, VectorizedArray< Number > > curl_in, const unsigned int q_point)
 
- Public Member Functions inherited from FEEvaluationBase< dim, dim, Number >
 ~FEEvaluationBase ()
 
void reinit (const unsigned int cell)
 
void reinit (const TriaIterator< DoFCellAccessor< DoFHandlerType, level_dof_access > > &cell)
 
void reinit (const typename Triangulation< dim >::cell_iterator &cell)
 
unsigned int get_cell_data_number () const
 
internal::MatrixFreeFunctions::CellType get_cell_type () const
 
const internal::MatrixFreeFunctions::ShapeInfo< Number > & get_shape_info () const
 
void fill_JxW_values (AlignedVector< VectorizedArray< Number > > &JxW_values) const
 
void read_dof_values (const VectorType &src, const unsigned int first_index=0)
 
void read_dof_values_plain (const VectorType &src, const unsigned int first_index=0)
 
void distribute_local_to_global (VectorType &dst, const unsigned int first_index=0) const
 
void set_dof_values (VectorType &dst, const unsigned int first_index=0) const
 
value_type get_dof_value (const unsigned int dof) const
 
void submit_dof_value (const value_type val_in, const unsigned int dof)
 
value_type get_value (const unsigned int q_point) const
 
void submit_value (const value_type val_in, const unsigned int q_point)
 
gradient_type get_gradient (const unsigned int q_point) const
 
void submit_gradient (const gradient_type grad_in, const unsigned int q_point)
 
Tensor< 1, n_components_, Tensor< 2, dim, VectorizedArray< Number > > > get_hessian (const unsigned int q_point) const
 
gradient_type get_hessian_diagonal (const unsigned int q_point) const
 
value_type get_laplacian (const unsigned int q_point) const
 
VectorizedArray< Number > get_divergence (const unsigned int q_point) const
 
SymmetricTensor< 2, dim, VectorizedArray< Number > > get_symmetric_gradient (const unsigned int q_point) const
 
Tensor< 1,(dim==2?1:dim), VectorizedArray< Number > > get_curl (const unsigned int q_point) const
 
void submit_divergence (const VectorizedArray< Number > div_in, const unsigned int q_point)
 
void submit_symmetric_gradient (const SymmetricTensor< 2, dim, VectorizedArray< Number > > grad_in, const unsigned int q_point)
 
void submit_curl (const Tensor< 1, dim==2?1:dim, VectorizedArray< Number > > curl_in, const unsigned int q_point)
 
value_type integrate_value () const
 
VectorizedArray< Number > JxW (const unsigned int q_point) const
 
const VectorizedArray< Number > * begin_dof_values () const
 
VectorizedArray< Number > * begin_dof_values ()
 
const VectorizedArray< Number > * begin_values () const
 
VectorizedArray< Number > * begin_values ()
 
const VectorizedArray< Number > * begin_gradients () const
 
VectorizedArray< Number > * begin_gradients ()
 
const VectorizedArray< Number > * begin_hessians () const
 
VectorizedArray< Number > * begin_hessians ()
 
const std::vector< unsigned int > & get_internal_dof_numbering () const
 
ArrayView< VectorizedArray< Number > > get_scratch_data () const
 

Protected Member Functions

 FEEvaluationAccess (const MatrixFree< dim, Number > &matrix_free, const unsigned int fe_no, const unsigned int quad_no, const unsigned int dofs_per_cell, const unsigned int n_q_points)
 
template<int n_components_other>
 FEEvaluationAccess (const Mapping< dim > &mapping, const FiniteElement< dim > &fe, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags, const unsigned int first_selected_component, const FEEvaluationBase< dim, n_components_other, Number > *other)
 
 FEEvaluationAccess (const FEEvaluationAccess &other)
 
FEEvaluationAccessoperator= (const FEEvaluationAccess &other)
 
- Protected Member Functions inherited from FEEvaluationBase< dim, dim, Number >
 FEEvaluationBase (const MatrixFree< dim, Number > &matrix_free, const unsigned int fe_no, const unsigned int quad_no, const unsigned int fe_degree, const unsigned int n_q_points)
 
 FEEvaluationBase (const Mapping< dim > &mapping, const FiniteElement< dim > &fe, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags, const unsigned int first_selected_component, const FEEvaluationBase< dim, n_components_other, Number > *other)
 
 FEEvaluationBase (const FEEvaluationBase &other)
 
void read_dof_values_plain (const VectorType *src_data[])
 
FEEvaluationBaseoperator= (const FEEvaluationBase &other)
 
void read_write_operation (const VectorOperation &operation, VectorType *vectors[]) const
 

Additional Inherited Members

- Protected Attributes inherited from FEEvaluationBase< dim, dim, Number >
AlignedVector< VectorizedArray< Number > > * scratch_data_array
 
VectorizedArray< Number > * scratch_data
 
VectorizedArray< Number > * values_dofs [n_components]
 
VectorizedArray< Number > * values_quad [n_components]
 
VectorizedArray< Number > * gradients_quad [n_components][dim]
 
VectorizedArray< Number > * hessians_quad [n_components][(dim *(dim+1))/2]
 
const unsigned int quad_no
 
const unsigned int n_fe_components
 
const unsigned int active_fe_index
 
const unsigned int active_quad_index
 
const MatrixFree< dim, Number > * matrix_info
 
const internal::MatrixFreeFunctions::DoFInfo * dof_info
 
const internal::MatrixFreeFunctions::MappingInfo< dim, Number > * mapping_info
 
const internal::MatrixFreeFunctions::ShapeInfo< Number > * data
 
const Tensor< 1, dim, VectorizedArray< Number > > * cartesian_data
 
const Tensor< 2, dim, VectorizedArray< Number > > * jacobian
 
const VectorizedArray< Number > * J_value
 
const VectorizedArray< Number > * quadrature_weights
 
const Point< dim, VectorizedArray< Number > > * quadrature_points
 
const Tensor< 2, dim, VectorizedArray< Number > > * jacobian_grad
 
const Tensor< 1,(dim >1?dim *(dim-1)/2:1), Tensor< 1, dim, VectorizedArray< Number > > > * jacobian_grad_upper
 
unsigned int cell
 
internal::MatrixFreeFunctions::CellType cell_type
 
unsigned int cell_data_number
 
bool dof_values_initialized
 
bool values_quad_initialized
 
bool gradients_quad_initialized
 
bool hessians_quad_initialized
 
bool values_quad_submitted
 
bool gradients_quad_submitted
 
std_cxx1x::shared_ptr< internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number > > mapped_geometry
 
std::vector< types::global_dof_indexold_style_dof_indices
 
const unsigned int first_selected_component
 
std::vector< types::global_dof_indexlocal_dof_indices
 

Detailed Description

template<int dim, typename Number>
class FEEvaluationAccess< dim, dim, Number >

This class provides access to the data fields of the FEEvaluation classes. Partial specialization for fields with as many components as the underlying space dimension, i.e., values are of type Tensor<1,dim> and gradients of type Tensor<2,dim>. Provides some additional functions for access, like the symmetric gradient and divergence.

Author
Katharina Kormann and Martin Kronbichler, 2010, 2011

Definition at line 1161 of file fe_evaluation.h.

Constructor & Destructor Documentation

◆ FEEvaluationAccess() [1/3]

template<int dim, typename Number >
FEEvaluationAccess< dim, dim, Number >::FEEvaluationAccess ( const MatrixFree< dim, Number > &  matrix_free,
const unsigned int  fe_no,
const unsigned int  quad_no,
const unsigned int  dofs_per_cell,
const unsigned int  n_q_points 
)
protected

Constructor. Made protected to avoid initialization in user code. Takes all data stored in MatrixFree. If applied to problems with more than one finite element or more than one quadrature formula selected during construction of matrix_free, fe_no and quad_no allow to select the appropriate components.

◆ FEEvaluationAccess() [2/3]

template<int dim, typename Number >
template<int n_components_other>
FEEvaluationAccess< dim, dim, Number >::FEEvaluationAccess ( const Mapping< dim > &  mapping,
const FiniteElement< dim > &  fe,
const Quadrature< 1 > &  quadrature,
const UpdateFlags  update_flags,
const unsigned int  first_selected_component,
const FEEvaluationBase< dim, n_components_other, Number > *  other 
)
protected

Constructor with reduced functionality for similar usage of FEEvaluation as FEValues, including matrix assembly.

◆ FEEvaluationAccess() [3/3]

template<int dim, typename Number >
FEEvaluationAccess< dim, dim, Number >::FEEvaluationAccess ( const FEEvaluationAccess< dim, dim, Number > &  other)
protected

Copy constructor

Member Function Documentation

◆ get_gradient()

template<int dim, typename Number >
gradient_type FEEvaluationAccess< dim, dim, Number >::get_gradient ( const unsigned int  q_point) const

Return the gradient of a finite element function at quadrature point number q_point after a call to evaluate(...,true,...).

◆ get_divergence()

template<int dim, typename Number >
VectorizedArray<Number> FEEvaluationAccess< dim, dim, Number >::get_divergence ( const unsigned int  q_point) const

Return the divergence of a vector-valued finite element at quadrature point number q_point after a call to evaluate(...,true,...).

◆ get_symmetric_gradient()

template<int dim, typename Number >
SymmetricTensor<2,dim,VectorizedArray<Number> > FEEvaluationAccess< dim, dim, Number >::get_symmetric_gradient ( const unsigned int  q_point) const

Return the symmetric gradient of a vector-valued finite element at quadrature point number q_point after a call to evaluate(...,true,...). It corresponds to 0.5 (grad+gradT).

◆ get_curl()

template<int dim, typename Number >
Tensor<1,(dim==2?1:dim),VectorizedArray<Number> > FEEvaluationAccess< dim, dim, Number >::get_curl ( const unsigned int  q_point) const

Return the curl of the vector field, \(\nabla \times v\) after a call to evaluate(...,true,...).

◆ get_hessian()

template<int dim, typename Number >
Tensor<3,dim,VectorizedArray<Number> > FEEvaluationAccess< dim, dim, Number >::get_hessian ( const unsigned int  q_point) const

Return the Hessian of a finite element function at quadrature point number q_point after a call to evaluate(...,true). If only the diagonal of the Hessian or its trace, the Laplacian, is needed, use the respective functions.

◆ get_hessian_diagonal()

template<int dim, typename Number >
gradient_type FEEvaluationAccess< dim, dim, Number >::get_hessian_diagonal ( const unsigned int  q_point) const

Return the diagonal of the Hessian of a finite element function at quadrature point number q_point after a call to evaluate(...,true).

◆ submit_gradient() [1/2]

template<int dim, typename Number >
void FEEvaluationAccess< dim, dim, Number >::submit_gradient ( const gradient_type  grad_in,
const unsigned int  q_point 
)

Write a contribution that is tested by the gradient to the field containing the values on quadrature points with component q_point. Access to the same field as through get_gradient. If applied before the function integrate(...,true) is called, this specifies what is tested by all basis function gradients on the current cell and integrated over.

◆ submit_gradient() [2/2]

template<int dim, typename Number >
void FEEvaluationAccess< dim, dim, Number >::submit_gradient ( const Tensor< 1, dim, Tensor< 1, dim, VectorizedArray< Number > > >  grad_in,
const unsigned int  q_point 
)

Write a contribution that is tested by the gradient to the field containing the values on quadrature points with component q_point. This function is an alternative to the other submit_gradient function when using a system of fixed number of equations which happens to coincide with the dimension for some dimensions, but not all. To allow for dimension-independent programming, this function can be used instead.

◆ submit_divergence()

template<int dim, typename Number >
void FEEvaluationAccess< dim, dim, Number >::submit_divergence ( const VectorizedArray< Number >  div_in,
const unsigned int  q_point 
)

Write a contribution that is tested by the divergence to the field containing the values on quadrature points with component q_point. Access to the same field as through get_gradient. If applied before the function integrate(...,true) is called, this specifies what is tested by all basis function gradients on the current cell and integrated over.

◆ submit_symmetric_gradient()

template<int dim, typename Number >
void FEEvaluationAccess< dim, dim, Number >::submit_symmetric_gradient ( const SymmetricTensor< 2, dim, VectorizedArray< Number > >  grad_in,
const unsigned int  q_point 
)

Write a contribution that is tested by the symmetric gradient to the field containing the values on quadrature points with component q_point. Access to the same field as through get_symmetric_gradient. If applied before the function integrate(...,true) is called, this specifies the symmetric gradient which is tested by all basis function symmetric gradients on the current cell and integrated over.

◆ submit_curl()

template<int dim, typename Number >
void FEEvaluationAccess< dim, dim, Number >::submit_curl ( const Tensor< 1, dim==2?1:dim, VectorizedArray< Number > >  curl_in,
const unsigned int  q_point 
)

Write the components of a curl containing the values on quadrature point q_point. Access to the same data field as through get_gradient.

◆ operator=()

template<int dim, typename Number >
FEEvaluationAccess& FEEvaluationAccess< dim, dim, Number >::operator= ( const FEEvaluationAccess< dim, dim, Number > &  other)
protected

Copy assignment operator


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