Reference documentation for deal.II version 9.1.1
|
#include <deal.II/matrix_free/cuda_fe_evaluation.h>
Public Member Functions | |
FEEvaluation (const unsigned int cell_id, 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 bool evaluate_val, const bool evaluate_grad) |
void | integrate (const bool integrate_val, const bool integrate_grad) |
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) |
template<typename Functor > | |
void | apply_quad_point_operations (const Functor &func) |
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 72 of file cuda_fe_evaluation.h.
FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::FEEvaluation | ( | const unsigned int | cell_id, |
const data_type * | data, | ||
SharedData< dim, Number > * | shdata | ||
) |
Constructor.
Definition at line 203 of file cuda_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_valuess as well.
Definition at line 228 of file cuda_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 257 of file cuda_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 282 of file cuda_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 316 of file cuda_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 | ( | const unsigned int | q_point | ) | const |
Return the value of a finite element function at quadrature point number q_point
after a call to evalue
(true,...).
Definition at line 357 of file cuda_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::submit_value | ( | const value_type & | val_in, |
const unsigned int | q_point | ||
) |
Write a value to the field containing the values on quadrature points with component q_point
. Access to the same fields as through get_value()
, This specifies the value which is tested by all basis function on the current cell and integrated over.
Definition at line 372 of file cuda_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 | ( | 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).
Definition at line 390 of file cuda_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, 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
Definition at line 418 of file cuda_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::apply_quad_point_operations | ( | const Functor & | func | ) |
Apply the functor func
on every quadrature point.
func
needs to define
Definition at line 442 of file cuda_fe_evaluation.h.