![]() |
Reference documentation for deal.II version 8.4.2
|
#include <deal.II/matrix_free/fe_evaluation.h>
Public Member Functions | |
| FEEvaluation (const MatrixFree< dim, Number > &matrix_free, const unsigned int fe_no=0, const unsigned int quad_no=0) | |
| FEEvaluation (const Mapping< dim > &mapping, const FiniteElement< dim > &fe, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags, const unsigned int first_selected_component=0) | |
| FEEvaluation (const FiniteElement< dim > &fe, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags, const unsigned int first_selected_component=0) | |
| template<int n_components_other> | |
| FEEvaluation (const FiniteElement< dim > &fe, const FEEvaluationBase< dim, n_components_other, Number > &other, const unsigned int first_selected_component=0) | |
| FEEvaluation (const FEEvaluation &other) | |
| void | evaluate (const bool evaluate_val, const bool evaluate_grad, const bool evaluate_hess=false) |
| void | integrate (const bool integrate_val, const bool integrate_grad) |
| Point< dim, VectorizedArray< Number > > | quadrature_point (const unsigned int q_point) const |
Public Member Functions inherited from FEEvaluationBase< dim, n_components_, Number > | |
| void | reinit (const unsigned int cell) |
| template<typename DoFHandlerType , bool level_dof_access> | |
| 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 |
| template<typename VectorType > | |
| void | read_dof_values (const VectorType &src) |
| template<typename VectorType > | |
| void | read_dof_values (const std::vector< VectorType > &src, const unsigned int first_index=0) |
| template<typename VectorType > | |
| void | read_dof_values (const std::vector< VectorType *> &src, const unsigned int first_index=0) |
| template<typename VectorType > | |
| void | read_dof_values_plain (const VectorType &src) |
| template<typename VectorType > | |
| void | read_dof_values_plain (const std::vector< VectorType > &src, const unsigned int first_index=0) |
| template<typename VectorType > | |
| void | read_dof_values_plain (const std::vector< VectorType *> &src, const unsigned int first_index=0) |
| template<typename VectorType > | |
| void | distribute_local_to_global (VectorType &dst) const |
| template<typename VectorType > | |
| void | distribute_local_to_global (std::vector< VectorType > &dst, const unsigned int first_index=0) const |
| template<typename VectorType > | |
| void | distribute_local_to_global (std::vector< VectorType *> &dst, const unsigned int first_index=0) const |
| template<typename VectorType > | |
| void | set_dof_values (VectorType &dst) const |
| template<typename VectorType > | |
| void | set_dof_values (std::vector< VectorType > &dst, const unsigned int first_index=0) const |
| template<typename VectorType > | |
| 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 () |
| const std::vector< unsigned int > & | get_internal_dof_numbering () const |
Public Attributes | |
| const unsigned int | dofs_per_cell |
Private Member Functions | |
| void | check_template_arguments (const unsigned int fe_no) |
| void | set_data_pointers () |
Private Attributes | |
| VectorizedArray< Number > | my_data_array [n_components *(tensor_dofs_per_cell+1+(dim *dim+2 *dim+1) *n_q_points)] |
| void(* | evaluate_funct )(const internal::MatrixFreeFunctions::ShapeInfo< Number > &, VectorizedArray< Number > *values_dofs_actual[], VectorizedArray< Number > *values_quad[], VectorizedArray< Number > *gradients_quad[][dim], VectorizedArray< Number > *hessians_quad[][(dim *(dim+1))/2], const bool evaluate_val, const bool evaluate_grad, const bool evaluate_lapl) |
| void(* | integrate_funct )(const internal::MatrixFreeFunctions::ShapeInfo< Number > &, VectorizedArray< Number > *values_dofs_actual[], VectorizedArray< Number > *values_quad[], VectorizedArray< Number > *gradients_quad[][dim], const bool evaluate_val, const bool evaluate_grad) |
Additional Inherited Members | |
Protected Member Functions inherited from FEEvaluationAccess< dim, n_components_, Number > | |
| FEEvaluationAccess (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) | |
| 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) | |
Protected Member Functions inherited from FEEvaluationBase< dim, n_components_, 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) | |
| 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) | |
| template<typename VectorType , typename VectorOperation > | |
| void | read_write_operation (const VectorOperation &operation, VectorType *vectors[]) const |
| template<typename VectorType > | |
| void | read_dof_values_plain (const VectorType *src_data[]) |
Protected Attributes inherited from FEEvaluationBase< dim, n_components_, Number > | |
| 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 |
| std_cxx11::shared_ptr< internal::MatrixFreeFunctions::ShapeInfo< Number > > | stored_shape_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_index > | old_style_dof_indices |
| const unsigned int | first_selected_component |
| std::vector< types::global_dof_index > | local_dof_indices |
The class that provides all functions necessary to evaluate functions at quadrature points and cell integrations. In functionality, this class is similar to FEValues<dim>, however, it includes a lot of specialized functions that make it much faster (between 5 and 500, depending on the polynomial order).
The first and foremost way of usage is to initialize this class from a MatrixFree object that caches everything related to the degrees of freedom and the mapping information. This way, it is possible to use vectorization for applying a vector operation for several cells at once. This setting is explained in the step-37 and step-48 tutorial programs. For vector-valued problems, the deal.II test suite includes a few additional examples as well, e.g. the Stokes operator found at https://github.com/dealii/dealii/blob/master/tests/matrix_free/matrix_vector_stokes_noflux.cc
For most operator evaluation tasks, this path provides the most efficient solution by combining pre-computed data for the mapping (Jacobian transformations for the geometry description) with on-the-fly evaluation of basis functions. In other words, the framework provides a trade-off between memory usage and initialization of objects that is suitable for matrix-free operator evaluation.
The second form of usage is to initialize FEEvaluation from geometry information generated by FEValues. This allows to apply the integration loops on the fly without prior initialization of MatrixFree objects. This can be useful when the memory and initialization cost of MatrixFree is not acceptable, e.g. when a different number of quadrature points should be used for one single evaluation in error computation. Also, when using the routines of this class to assemble matrices the trade-off implied by the MatrixFree class may not be desired. In such a case, the cost to initialize the necessary geometry data on the fly is comparably low and thus avoiding a global object MatrixFree can be useful. When used in this way, reinit methods reminiscent from FEValues with a cell iterator are to be used. However, note that this model results in working on a single cell at a time, with geometry data duplicated in all components of the vectorized array. Thus, vectorization is only useful when it can apply the same operation on different data, e.g. when performing matrix assembly.
As an example, consider the following code to assemble the contributions to the Laplace matrix:
This code generates the columns of the cell matrix with the loop over i above. The way this is done is the following: FEEvaluation's routines focus on the evaluation of finite element operators, so the way to get a cell matrix out of an operator evaluation is to apply it to all the unit vectors on the cell. This might seem inefficient but the evaluation routines used here are so quick that they still work much faster than what is possible with FEValues.
Due to vectorization, we can actually generate matrix data for several unit vectors at a time (e.g. 4). The variable n_items make sure that we do the last iteration where the number of cell dofs is not divisible by the vectorization length correctly. Also note that we need to get the internal dof numbering applied by fe_eval because FEEvaluation internally uses a lexicographic numbering of degrees of freedom. This is necessary to efficiently work with tensor products where all degrees of freedom along a dimension must be laid out in a regular way.
This class contains specialized evaluation routines for several elements based on tensor-product quadrature formulas and tensor-product-like shape functions, including standard FE_Q or FE_DGQ elements and quadrature points symmetric around 0.5 (like Gauss quadrature), FE_DGP elements based on truncated tensor products as well as the faster case of Gauss-Lobatto elements with Gauss-Lobatto quadrature which give diagonal mass matrices and quicker evaluation internally. The main benefit of this class is the evaluation of all shape functions in all quadrature or integration over all shape functions in dim (fe_degree+1)dim+1 operations instead of the slower (fe_degree+1)2*dim complexity in the evaluation routines of FEValues.
Note that many of the operations available through this class are inherited from the base class FEEvaluationBase, in particular reading from and writing to vectors. Also, the class inherits from FEEvaluationAccess that implements access to values, gradients and Hessians of the finite element function on quadrature points.
This class assumes that shape functions of the FiniteElement under consideration do not depend on the actual shape of the cells in real space. Currently, other finite elements cannot be treated with the matrix-free concept.
This class has five template arguments:
| 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 |
| n_q_points_1d | Number of points in the quadrature formula in 1D, defaults to 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). Defaults to 1. |
| Number | Number format, usually double or float. Defaults to double |
Definition at line 50 of file fe_evaluation.h.
| FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::FEEvaluation | ( | const MatrixFree< dim, Number > & | matrix_free, |
| const unsigned int | fe_no = 0, |
||
| const unsigned int | quad_no = 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, fe_no and quad_no allow to select the appropriate components.
| FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::FEEvaluation | ( | const Mapping< dim > & | mapping, |
| const FiniteElement< dim > & | fe, | ||
| const Quadrature< 1 > & | quadrature, | ||
| const UpdateFlags | update_flags, | ||
| const unsigned int | first_selected_component = 0 |
||
| ) |
Constructor that comes with reduced functionality and works similar as FEValues. The arguments are similar to the ones passed to the constructor of FEValues, with the notable difference that FEEvaluation expects a one- dimensional quadrature formula, Quadrature<1>, instead of a dim dimensional one. The finite element can be both scalar or vector valued, but this method always only selects a scalar base element at a time (with n_components copies as specified by the class template). For vector- valued elements, the optional argument first_selected_component allows to specify the index of the base element to be used for evaluation. Note that the internal data structures always assume that the base element is primitive, non-primitive are not supported currently.
As known from FEValues, a call to the reinit method with a Triangulation<dim>::cell_iterator is necessary to make the geometry and degrees of freedom of the current class known. If the iterator includes DoFHandler information (i.e., it is a DoFHandler<dim>::cell_iterator or similar), the initialization allows to also read from or write to vectors in the standard way for DoFHandler<dim>::active_cell_iterator types for one cell at a time. However, this approach is much slower than the path with MatrixFree with MPI since index translation has to be done. As only one cell at a time is used, this method does not vectorize over several elements (which is most efficient for vector operations), but only possibly within the element if the evaluate/integrate routines are combined inside user code (e.g. for computing cell matrices).
| FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::FEEvaluation | ( | const FiniteElement< dim > & | fe, |
| const Quadrature< 1 > & | quadrature, | ||
| const UpdateFlags | update_flags, | ||
| const unsigned int | first_selected_component = 0 |
||
| ) |
Constructor for the reduced functionality. This constructor is equivalent to the other one except that it makes the object use a
mapping (i.e., an object of type MappingQGeneric(1)) implicitly.
| FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::FEEvaluation | ( | const FiniteElement< dim > & | fe, |
| const FEEvaluationBase< dim, n_components_other, Number > & | other, | ||
| const unsigned int | first_selected_component = 0 |
||
| ) |
Constructor for the reduced functionality. Similar to the other constructor with FiniteElement argument but using another FEEvaluationBase object to provide info about the geometry. This allows several FEEvaluation objects to share the geometry evaluation, i.e., the underlying mapping and quadrature points do only need to be evaluated once. Make sure to not pass an optional object around when you intend to use the FEEvaluation object in parallel to the given one because otherwise the intended sharing may create race conditions.
| FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::FEEvaluation | ( | const FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number > & | other | ) |
Copy constructor. If FEEvaluationBase was constructed from a mapping, fe, quadrature, and update flags, the underlying geometry evaluation based on FEValues will be deep-copied in order to allow for using in parallel with threads.
| void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::evaluate | ( | const bool | evaluate_val, |
| const bool | evaluate_grad, | ||
| const bool | evaluate_hess = false |
||
| ) |
Evaluates the function values, the gradients, and the Laplacians of the FE function given at the DoF values in the input vector 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_laplacian give useful information (unless these values have been set manually).
| void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::integrate | ( | const bool | integrate_val, |
| const bool | integrate_grad | ||
| ) |
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.
| Point<dim,VectorizedArray<Number> > FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::quadrature_point | ( | const unsigned int | q_point | ) | const |
Returns the q-th quadrature point stored in MappingInfo.
|
private |
Checks if the template arguments regarding degree of the element corresponds to the actual element used at initialization.
|
private |
Sets the pointers of the base class to my_data_array.
| const unsigned int FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::dofs_per_cell |
The number of scalar degrees of freedom on the cell. Usually close to tensor_dofs_per_cell, but depends on the actual element selected and is thus not static.
Definition at line 1650 of file fe_evaluation.h.
|
private |
Internally stored variables for the different data fields.
Definition at line 1656 of file fe_evaluation.h.
|
private |
Function pointer for the evaluate function
Definition at line 1672 of file fe_evaluation.h.
|
private |
Function pointer for the integrate function
Definition at line 1684 of file fe_evaluation.h.
1.8.12