![]() |
Reference documentation for deal.II version 8.5.1
|
#include <deal.II/matrix_free/fe_evaluation.h>
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) |
![]() | |
~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) | |
FEEvaluationAccess & | operator= (const FEEvaluationAccess &other) |
![]() | |
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[]) |
FEEvaluationBase & | operator= (const FEEvaluationBase &other) |
void | read_write_operation (const VectorOperation &operation, VectorType *vectors[]) const |
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.
Definition at line 1161 of file fe_evaluation.h.
|
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.
|
protected |
Constructor with reduced functionality for similar usage of FEEvaluation as FEValues, including matrix assembly.
|
protected |
Copy constructor
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,...).
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,...).
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)
.
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,...).
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.
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).
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.
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.
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.
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.
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
.
|
protected |
Copy assignment operator