Reference documentation for deal.II version 8.1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
Public Types | Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
FEEvaluationAccess< dim, dofs_per_cell_, n_q_points_, dim, Number > Class Template Reference

#include <fe_evaluation.h>

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

Public Types

typedef Number number_type
 
typedef Tensor< 1, dim,
VectorizedArray< Number > > 
value_type
 
typedef Tensor< 2, dim,
VectorizedArray< Number > > 
gradient_type
 
typedef FEEvaluationBase< dim,
dofs_per_cell_, n_q_points_,
dim, Number > 
BaseClass
 
- Public Types inherited from FEEvaluationBase< dim, dofs_per_cell_, n_q_points_, dim, Number >
typedef Number number_type
 
typedef Tensor
< 1, n_components_,
VectorizedArray< Number > > 
value_type
 
typedef Tensor
< 1, n_components_, Tensor
< 1, dim, VectorizedArray
< Number > > > 
gradient_type
 

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, dofs_per_cell_, n_q_points_, dim, Number >
void reinit (const unsigned int cell)
 
unsigned int get_cell_data_number () const
 
internal::MatrixFreeFunctions::CellType get_cell_type () const
 
void read_dof_values (const VectorType &src)
 
void read_dof_values (const std::vector< VectorType > &src, const unsigned int first_index=0)
 
void read_dof_values (const std::vector< VectorType * > &src, const unsigned int first_index=0)
 
void read_dof_values_plain (const VectorType &src)
 
void read_dof_values_plain (const std::vector< VectorType > &src, const unsigned int first_index=0)
 
void read_dof_values_plain (const std::vector< VectorType * > &src, const unsigned int first_index=0)
 
void distribute_local_to_global (VectorType &dst) const
 
void distribute_local_to_global (std::vector< VectorType > &dst, const unsigned int first_index=0) const
 
void distribute_local_to_global (std::vector< VectorType * > &dst, const unsigned int first_index=0) const
 
void set_dof_values (VectorType &dst) const
 
void set_dof_values (std::vector< VectorType > &dst, const unsigned int first_index=0) const
 
void set_dof_values (std::vector< 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
 
value_type integrate_value () 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 ()
 

Static Public Attributes

static const unsigned int dimension = dim
 
static const unsigned int n_components = dim
 
static const unsigned int dofs_per_cell = dofs_per_cell_
 
static const unsigned int n_q_points = n_q_points_
 
- Static Public Attributes inherited from FEEvaluationBase< dim, dofs_per_cell_, n_q_points_, dim, Number >
static const unsigned int dimension
 
static const unsigned int n_components
 
static const unsigned int dofs_per_cell
 
static const unsigned int n_q_points
 

Protected Member Functions

 FEEvaluationAccess (const MatrixFree< dim, Number > &matrix_free, const unsigned int fe_no=0, const unsigned int quad_no=0)
 
- Protected Member Functions inherited from FEEvaluationBase< dim, dofs_per_cell_, n_q_points_, dim, Number >
void read_dof_values_plain (const VectorType *src_data[])
 
 FEEvaluationBase (const MatrixFree< dim, Number > &matrix_free, const unsigned int fe_no=0, const unsigned int quad_no=0)
 
void read_write_operation (const VectorOperation &operation, VectorType *vectors[]) const
 

Additional Inherited Members

- Protected Attributes inherited from FEEvaluationBase< dim, dofs_per_cell_, n_q_points_, dim, Number >
VectorizedArray< Number > values_dofs [n_components][dofs_per_cell >0?dofs_per_cell:1]
 
VectorizedArray< Number > values_quad [n_components][n_q_points >0?n_q_points:1]
 
VectorizedArray< Number > gradients_quad [n_components][dim][n_q_points >0?n_q_points:1]
 
VectorizedArray< Number > hessians_quad [n_components][(dim *(dim+1))/2][n_q_points >0?n_q_points:1]
 
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
 

Detailed Description

template<int dim, int dofs_per_cell_, int n_q_points_, typename Number>
class FEEvaluationAccess< dim, dofs_per_cell_, n_q_points_, 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 911 of file fe_evaluation.h.

Constructor & Destructor Documentation

template<int dim, int dofs_per_cell_, int n_q_points_, typename Number >
FEEvaluationAccess< dim, dofs_per_cell_, n_q_points_, dim, Number >::FEEvaluationAccess ( const MatrixFree< dim, Number > &  matrix_free,
const unsigned int  fe_no = 0,
const unsigned int  quad_no = 0 
)
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.

Member Function Documentation

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

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

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

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

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

Returns 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).

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

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

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

Returns 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.

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

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

template<int dim, int dofs_per_cell_, int n_q_points_, typename Number >
void FEEvaluationAccess< dim, dofs_per_cell_, n_q_points_, 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.

template<int dim, int dofs_per_cell_, int n_q_points_, typename Number >
void FEEvaluationAccess< dim, dofs_per_cell_, n_q_points_, 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.

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

Write a constribution 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.

template<int dim, int dofs_per_cell_, int n_q_points_, typename Number >
void FEEvaluationAccess< dim, dofs_per_cell_, n_q_points_, 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 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 the gradient which is tested by all basis function gradients on the current cell and integrated over.

template<int dim, int dofs_per_cell_, int n_q_points_, typename Number >
void FEEvaluationAccess< dim, dofs_per_cell_, n_q_points_, 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.


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