Reference documentation for deal.II version 9.0.0
|
#include <deal.II/matrix_free/fe_evaluation.h>
Public Types | |
typedef FEEvaluationAccess< dim, n_components_, Number, true > | BaseClass |
typedef Number | number_type |
typedef BaseClass::value_type | value_type |
typedef BaseClass::gradient_type | gradient_type |
Public Member Functions | |
FEFaceEvaluation (const MatrixFree< dim, Number > &matrix_free, const bool is_interior_face=true, const unsigned int dof_no=0, const unsigned int quad_no=0, const unsigned int first_selected_component=0) | |
~FEFaceEvaluation () | |
void | reinit (const unsigned int face_batch_number) |
void | reinit (const unsigned int cell_batch_number, const unsigned int face_number) |
void | evaluate (const bool evaluate_values, const bool evaluate_gradients) |
void | evaluate (const VectorizedArray< Number > *values_array, const bool evaluate_values, const bool evaluate_gradients) |
template<typename VectorType > | |
void | gather_evaluate (const VectorType &input_vector, const bool evaluate_values, const bool evaluate_gradients) |
void | integrate (const bool integrate_values, const bool integrate_gradients) |
void | integrate (const bool integrate_values, const bool integrate_gradients, VectorizedArray< Number > *values_array) |
template<typename VectorType > | |
void | integrate_scatter (const bool integrate_values, const bool integrate_gradients, VectorType &output_vector) |
Point< dim, VectorizedArray< Number > > | quadrature_point (const unsigned int q_point) const |
Public Member Functions inherited from FEEvaluationBase< dim, n_components_, Number, is_face > | |
~FEEvaluationBase () | |
unsigned int | get_cell_data_number () const |
unsigned int | get_mapping_data_index_offset () const |
internal::MatrixFreeFunctions::GeometryType | get_cell_type () const |
const internal::MatrixFreeFunctions::ShapeInfo< VectorizedArray< Number > > & | get_shape_info () const |
template<typename VectorType > | |
void | read_dof_values (const VectorType &src, const unsigned int first_index=0) |
template<typename VectorType > | |
void | read_dof_values_plain (const VectorType &src, const unsigned int first_index=0) |
template<typename VectorType > | |
void | distribute_local_to_global (VectorType &dst, const unsigned int first_index=0) const |
template<typename VectorType > | |
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 |
value_type | get_normal_derivative (const unsigned int q_point) const |
void | submit_gradient (const gradient_type grad_in, const unsigned int q_point) |
void | submit_normal_derivative (const value_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_index) const |
void | fill_JxW_values (AlignedVector< VectorizedArray< Number > > &JxW_values) const |
Tensor< 2, dim, VectorizedArray< Number > > | inverse_jacobian (const unsigned int q_index) const |
Tensor< 1, dim, VectorizedArray< Number > > | get_normal_vector (const unsigned int q_point) const |
VectorizedArray< Number > | read_cell_data (const AlignedVector< VectorizedArray< Number > > &array) 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 |
Public Attributes | |
const unsigned int | dofs_per_component |
const unsigned int | dofs_per_cell |
const unsigned int | n_q_points |
Static Public Attributes | |
static constexpr unsigned int | dimension = dim |
static constexpr unsigned int | n_components = n_components_ |
static constexpr unsigned int | static_n_q_points = Utilities::pow(n_q_points_1d, dim-1) |
static constexpr unsigned int | static_n_q_points_cell = Utilities::pow(n_q_points_1d, dim) |
static constexpr unsigned int | static_dofs_per_component = Utilities::pow(fe_degree + 1, dim) |
static constexpr unsigned int | tensor_dofs_per_cell = static_dofs_per_component *n_components |
static constexpr unsigned int | static_dofs_per_cell = static_dofs_per_component *n_components |
Protected Member Functions | |
void | adjust_for_face_orientation (const bool integrate, const bool values, const bool gradients) |
Protected Member Functions inherited from FEEvaluationAccess< dim, n_components_, Number, true > | |
FEEvaluationAccess (const MatrixFree< dim, Number > &matrix_free, const unsigned int dof_no, const unsigned int first_selected_component, const unsigned int quad_no, const unsigned int fe_degree, const unsigned int n_q_points, const bool is_interior_face=true) | |
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, is_face > *other) | |
FEEvaluationAccess (const FEEvaluationAccess &other) | |
FEEvaluationAccess & | operator= (const FEEvaluationAccess &other) |
Protected Member Functions inherited from FEEvaluationBase< dim, n_components_, Number, is_face > | |
FEEvaluationBase (const MatrixFree< dim, Number > &matrix_free, const unsigned int dof_no, const unsigned int first_selected_component, const unsigned int quad_no, const unsigned int fe_degree, const unsigned int n_q_points, const bool is_interior_face) | |
template<int n_components_other> | |
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) | |
FEEvaluationBase & | operator= (const FEEvaluationBase &other) |
template<typename VectorType , typename VectorOperation > | |
void | read_write_operation (const VectorOperation &operation, VectorType *vectors[], const bool apply_constraints=true) const |
template<typename VectorType , typename VectorOperation > | |
void | read_write_operation_contiguous (const VectorOperation &operation, VectorType *vectors[]) const |
template<typename VectorType , typename VectorOperation > | |
void | read_write_operation_global (const VectorOperation &operation, VectorType *vectors[]) const |
The class that provides all functions necessary to evaluate functions at quadrature points and face integrations. The design of the class is similar to FEEvaluation and most of the interfaces are shared with that class, in particular most access functions that come from the common base classes FEEvaluationAccess and FEEvaluatioBase. Furthermore, the relation of this class to FEEvaluation is similar to the relation between FEValues and FEFaceValues.
dim | Dimension in which this class is to be used |
fe_degree | Degree of the tensor product finite element with fe_degree+1 degrees of freedom per coordinate direction. If set to -1, the degree of the underlying element will be used, which acts as a run time constant rather than a compile time constant that slows down the execution. |
n_q_points_1d | Number of points in the quadrature formula in 1D, usually chosen as fe_degree+1 |
n_components | Number of vector components when solving a system of PDEs. If the same operation is applied to several components of a PDE (e.g. a vector Laplace equation), they can be applied simultaneously with one call (and often more efficiently) |
Number | Number format, usually double or float |
Definition at line 2434 of file fe_evaluation.h.
typedef FEEvaluationAccess<dim,n_components_,Number,true> FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::BaseClass |
A typedef to the base class.
Definition at line 2440 of file fe_evaluation.h.
typedef Number FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::number_type |
A underlying number type specified as template argument.
Definition at line 2445 of file fe_evaluation.h.
typedef BaseClass::value_type FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::value_type |
The type of function values, e.g. VectorizedArray<Number>
for n_components=1
or Tensor<1,dim,VectorizedArray<Number> >
for n_components=dim
.
Definition at line 2452 of file fe_evaluation.h.
typedef BaseClass::gradient_type FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::gradient_type |
The type of gradients, e.g. Tensor<1,dim,VectorizedArray<Number>>
for n_components=1
or Tensor<2,dim,VectorizedArray<Number> >
for n_components=dim
.
Definition at line 2459 of file fe_evaluation.h.
FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::FEFaceEvaluation | ( | const MatrixFree< dim, Number > & | matrix_free, |
const bool | is_interior_face = true , |
||
const unsigned int | dof_no = 0 , |
||
const unsigned int | quad_no = 0 , |
||
const unsigned int | first_selected_component = 0 |
||
) |
Constructor. 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
, the appropriate component can be selected by the optional arguments.
matrix_free | Data object that contains all data |
is_interior_face | This selects which of the two cells of an internal face the current evaluator will be based upon. The interior face is the main face along which the normal vectors are oriented. The exterior face coming from the other side provides the same normal vector as the interior side, so if the outer normal vector to that side is desired, it must be multiplied by -1. |
dof_no | If matrix_free was set up with multiple DoFHandler objects, this parameter selects to which DoFHandler/ConstraintMatrix pair the given evaluator should be attached to. |
quad_no | If matrix_free was set up with multiple Quadrature objects, this parameter selects the appropriate number of the quadrature formula. |
first_selected_component | If the dof_handler selected by dof_no uses an FESystem consisting of more than one base element, this parameter selects the number of the base element in FESystem. Note that this does not directly relate to the component of the respective element due to the possibility for a multiplicity in the element. |
FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::~FEFaceEvaluation | ( | ) |
Destructor.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::reinit | ( | const unsigned int | face_batch_number | ) |
Initializes the operation pointer to the current face. This method is the default choice for face integration as the data stored in MappingInfo is stored according to this numbering. Unlike the reinit functions taking a cell iterator as argument below and the FEValues::reinit() methods, where the information related to a particular cell is generated in the reinit call, this function is very cheap since all data is pre-computed in matrix_free
, and only a few indices and pointers have to be set appropriately.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::reinit | ( | const unsigned int | cell_batch_number, |
const unsigned int | face_number | ||
) |
As opposed to the reinit() method from the base class, this reinit() method initializes for a given number of cells and a face number. This method is less efficient than the other reinit() method taking a numbering of the faces because it needs to copy the data associated with the faces to the cells in this call.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::evaluate | ( | const bool | evaluate_values, |
const bool | evaluate_gradients | ||
) |
Evaluates the function values, the gradients, and the Laplacians of the FE function given at the DoF values stored in the internal data field dof_values
(that is usually filled by the read_dof_values() method) at the quadrature points on the unit cell. The function arguments specify which parts shall actually be computed. Needs to be called before the functions get_value(), get_gradient() or get_normal_derivative() give useful information (unless these values have been set manually by accessing the internal data pointers).
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::evaluate | ( | const VectorizedArray< Number > * | values_array, |
const bool | evaluate_values, | ||
const bool | evaluate_gradients | ||
) |
Evaluates the function values, the gradients, and the Laplacians of the FE function given at the DoF values in the input array values_array
at the quadrature points on the unit cell. If multiple components are involved in the current FEEvaluation object, the sorting in values_array is such that all degrees of freedom for the first component come first, then all degrees of freedom for the second, and so on. The function arguments specify which parts shall actually be computed. Needs to be called before the functions get_value(), get_gradient(), or get_normal_derivative() give useful information (unless these values have been set manually).
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::gather_evaluate | ( | const VectorType & | input_vector, |
const bool | evaluate_values, | ||
const bool | evaluate_gradients | ||
) |
Reads from the input vector and evaluates the function values, the gradients, and the Laplacians of the FE function at the quadrature points on the unit cell. The function arguments specify which parts shall actually be computed. Needs to be called before the functions get_value(), get_gradient(), or get_normal_derivative() give useful information.
This call is equivalent to calling read_dof_values() followed by evaluate(), but might internally use some additional optimizations.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::integrate | ( | const bool | integrate_values, |
const bool | integrate_gradients | ||
) |
This function takes the values and/or gradients that are stored on quadrature points, tests them by all the basis functions/gradients on the cell and performs the cell integration. The two function arguments integrate_val
and integrate_grad
are used to enable/disable some of values or gradients. The result is written into the internal data field dof_values
(that is usually written into the result vector by the distribute_local_to_global() or set_dof_values() methods).
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::integrate | ( | const bool | integrate_values, |
const bool | integrate_gradients, | ||
VectorizedArray< Number > * | values_array | ||
) |
This function takes the values and/or gradients that are stored on quadrature points, tests them by all the basis functions/gradients on the cell and performs the cell integration. The two function arguments integrate_val
and integrate_grad
are used to enable/disable some of values or gradients. As opposed to the other integrate() method, this call stores the result of the testing in the given array values_array
.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::integrate_scatter | ( | const bool | integrate_values, |
const bool | integrate_gradients, | ||
VectorType & | output_vector | ||
) |
This function takes the values and/or gradients that are stored on quadrature points, tests them by all the basis functions/gradients on the cell and performs the cell integration. The two function arguments integrate_val
and integrate_grad
are used to enable/disable some of values or gradients.
This call is equivalent to calling integrate() followed by distribute_local_to_global(), but might internally use some additional optimizations.
Point<dim,VectorizedArray<Number> > FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::quadrature_point | ( | const unsigned int | q_point | ) | const |
Returns the q-th quadrature point on the face in real coordinates stored in MappingInfo.
|
protected |
For faces not oriented in the standard way, this method applies re-indexing on quadrature points. Called at the end of evaluate() and at the beginning of integrate().
|
static |
The dimension given as template argument.
Definition at line 2464 of file fe_evaluation.h.
|
static |
The number of solution components of the evaluator given as template argument.
Definition at line 2470 of file fe_evaluation.h.
|
static |
The static number of quadrature points determined from the given template argument n_q_points_1d
taken to the power of dim-1. Note that the actual number of quadrature points, n_q_points
, can be different if fe_degree=-1
is given and run-time loop lengths are used rather than compile time ones.
Definition at line 2478 of file fe_evaluation.h.
|
static |
The static number of quadrature points on a cell with the same quadrature formula. Note that this value is only present for simpler comparison with the cell quadrature, as the actual number of points is given to a face by the static_n_q_points
variable.
Definition at line 2486 of file fe_evaluation.h.
|
static |
The static number of degrees of freedom of a scalar component determined from the given template argument fe_degree
. Note that the actual number of degrees of freedom dofs_per_component
can be different if fe_degree=-1
is given.
Definition at line 2494 of file fe_evaluation.h.
|
static |
The static number of degrees of freedom of all components determined from the given template argument fe_degree
. Note that the actual number of degrees of freedom dofs_per_cell
can be different if fe_degree=-1
is given.
Definition at line 2502 of file fe_evaluation.h.
|
static |
The static number of degrees of freedom of all components determined from the given template argument fe_degree
. Note that the actual number of degrees of freedom dofs_per_cell
can be different if fe_degree=-1
is given.
Definition at line 2510 of file fe_evaluation.h.
const unsigned int FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::dofs_per_component |
The number of degrees of freedom of a single component on the cell for the underlying evaluation object. Usually close to static_dofs_per_component, but the number depends on the actual element selected and is thus not static.
Definition at line 2672 of file fe_evaluation.h.
const unsigned int FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::dofs_per_cell |
The number of degrees of freedom on the cell accumulated over all components in the current evaluation object. Usually close to static_dofs_per_cell = static_dofs_per_component*n_components, but the number depends on the actual element selected and is thus not static.
Definition at line 2680 of file fe_evaluation.h.
const unsigned int FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::n_q_points |
The number of quadrature points in use. If the number of quadrature points in 1d is given as a template, this number is simply the dim-1
-th power of that value. If the element degree is set to -1 (dynamic selection of element degree), the static value of quadrature points is inaccurate and this value must be used instead.
Definition at line 2689 of file fe_evaluation.h.