deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
|
#include <deal.II/matrix_free/portable_fe_evaluation.h>
Public Types | |
using | value_type = std::conditional_t<(n_components_==1), Number, Tensor< 1, n_components_, Number > > |
using | gradient_type = std::conditional_t< n_components_==1, Tensor< 1, dim, Number >, std::conditional_t< n_components_==dim, Tensor< 2, dim, Number >, Tensor< 1, n_components_, Tensor< 1, dim, Number > > > > |
using | data_type = typename MatrixFree< dim, Number >::Data |
Public Member Functions | |
FEEvaluation (const data_type *data, SharedData< dim, Number > *shdata) | |
void | read_dof_values (const Number *src) |
void | distribute_local_to_global (Number *dst) const |
void | evaluate (const EvaluationFlags::EvaluationFlags evaluate_flag) |
void | evaluate (const bool evaluate_val, const bool evaluate_grad) |
void | integrate (const EvaluationFlags::EvaluationFlags integration_flag) |
void | integrate (const bool integrate_val, const bool integrate_grad) |
value_type | get_value (int q_point) const |
value_type | get_dof_value (int q_point) const |
void | submit_value (const value_type &val_in, int q_point) |
void | submit_dof_value (const value_type &val_in, int q_point) |
gradient_type | get_gradient (int q_point) const |
void | submit_gradient (const gradient_type &grad_in, int q_point) |
template<typename Functor > | |
void | apply_for_each_quad_point (const Functor &func) |
Static Public Attributes | |
static constexpr unsigned int | dimension = dim |
static constexpr unsigned int | n_components = n_components_ |
static constexpr unsigned int | n_q_points |
static constexpr unsigned int | tensor_dofs_per_cell |
Private Attributes | |
const data_type * | data |
SharedData< dim, Number > * | shared_data |
int | cell_id |
This class provides all the functions necessary to evaluate functions at quadrature points and cell integrations. In functionality, this class is similar to FEValues<dim>.
This class has five template arguments:
dim | Dimension in which this class is to be used |
fe_degree | Degree of the tensor prodict finite element with fe_degree+1 degrees of freedom per coordinate direction |
n_q_points_1d | Number of points in the quadrature formular 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, double or float . Defaults to double . |
Definition at line 67 of file portable_fe_evaluation.h.
using Portable::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::value_type = std::conditional_t<(n_components_ == 1), Number, Tensor<1, n_components_, Number> > |
An alias for scalar quantities.
Definition at line 73 of file portable_fe_evaluation.h.
using Portable::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::gradient_type = std::conditional_t< n_components_ == 1, Tensor<1, dim, Number>, std::conditional_t<n_components_ == dim, Tensor<2, dim, Number>, Tensor<1, n_components_, Tensor<1, dim, Number> >> > |
An alias for vectorial quantities.
Definition at line 80 of file portable_fe_evaluation.h.
using Portable::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::data_type = typename MatrixFree<dim, Number>::Data |
An alias to kernel specific information.
Definition at line 90 of file portable_fe_evaluation.h.
FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::FEEvaluation | ( | const data_type * | data, |
SharedData< dim, Number > * | shdata | ||
) |
Constructor.
Definition at line 255 of file portable_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::read_dof_values | ( | const Number * | src | ) |
For the vector src
, read out the values on the degrees of freedom of the current cell, and store them internally. Similar functionality as the function DoFAccessor::get_interpolated_dof_values when no constraints are present, but it also includes constraints from hanging nodes, so once can see it as a similar function to AffineConstraints::read_dof_values() as well.
Definition at line 270 of file portable_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::distribute_local_to_global | ( | Number * | dst | ) | const |
Take the value stored internally on dof values of the current cell and sum them into the vector dst
. The function also applies constraints during the write operation. The functionality is hence similar to the function AffineConstraints::distribute_local_to_global.
Definition at line 301 of file portable_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::evaluate | ( | const EvaluationFlags::EvaluationFlags | evaluate_flag | ) |
Evaluate the function values and the gradients of the FE function given at the DoF values in the input vector at the quadrature points on the unit cell. The function argument evaluate_flag
specifies which parts shall actually be computed. This function needs to be called before the functions get_value()
or get_gradient()
give useful information.
Definition at line 345 of file portable_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::evaluate | ( | const bool | evaluate_val, |
const bool | evaluate_grad | ||
) |
Evaluate the function values and the gradients 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. This function needs to be called before the functions get_value()
or get_gradient()
give useful information.
Definition at line 397 of file portable_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::integrate | ( | const EvaluationFlags::EvaluationFlags | integration_flag | ) |
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 as specified by the integration_flag
argument.
Definition at line 414 of file portable_fe_evaluation.h.
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 the values or the gradients.
Definition at line 464 of file portable_fe_evaluation.h.
FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::value_type FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::get_value | ( | int | q_point | ) | const |
Same as above, except that the quadrature point is computed from thread id.
Definition at line 485 of file portable_fe_evaluation.h.
FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::value_type FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::get_dof_value | ( | int | q_point | ) | const |
Same as above, except that the local dof index is computed from the thread id.
Definition at line 513 of file portable_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::submit_value | ( | const value_type & | val_in, |
int | q_point | ||
) |
Same as above, except that the quadrature point is computed from the thread id.
Definition at line 537 of file portable_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::submit_dof_value | ( | const value_type & | val_in, |
int | q_point | ||
) |
Same as above, except that the local dof index is computed from the thread id.
Definition at line 560 of file portable_fe_evaluation.h.
FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::gradient_type FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::get_gradient | ( | int | q_point | ) | const |
Same as above, except that the quadrature point is computed from the thread id.
Definition at line 586 of file portable_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::submit_gradient | ( | const gradient_type & | grad_in, |
int | q_point | ||
) |
Same as above, except that the quadrature point is computed from the thread id.
Definition at line 626 of file portable_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::apply_for_each_quad_point | ( | const Functor & | func | ) |
Same as above, except that the functor func
only takes a single input argument (fe_eval) and computes the quadrature point from the thread id.
func
needs to define
Definition at line 665 of file portable_fe_evaluation.h.
|
staticconstexpr |
Dimension.
Definition at line 95 of file portable_fe_evaluation.h.
|
staticconstexpr |
Number of components.
Definition at line 100 of file portable_fe_evaluation.h.
|
staticconstexpr |
Number of quadrature points per cell.
Definition at line 105 of file portable_fe_evaluation.h.
|
staticconstexpr |
Number of tensor degrees of freedoms per cell.
Definition at line 111 of file portable_fe_evaluation.h.
|
private |
Definition at line 242 of file portable_fe_evaluation.h.
|
private |
Definition at line 243 of file portable_fe_evaluation.h.
|
private |
Definition at line 244 of file portable_fe_evaluation.h.