Reference documentation for deal.II version 9.3.3
|
#include <deal.II/matrix_free/cuda_fe_evaluation.h>
Public Types | |
using | value_type = Number |
using | gradient_type = Tensor< 1, dim, Number > |
using | data_type = typename MatrixFree< dim, Number >::Data |
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 |
value_type | get_value () const |
value_type | get_dof_value (const unsigned int dof) const |
value_type | get_dof_value () const |
void | submit_value (const value_type &val_in, const unsigned int q_point) |
void | submit_value (const value_type &val_in) |
void | submit_dof_value (const value_type &val_in, const unsigned int dof) |
void | submit_dof_value (const value_type &val_in) |
gradient_type | get_gradient (const unsigned int q_point) const |
gradient_type | get_gradient () const |
void | submit_gradient (const gradient_type &grad_in, const unsigned int q_point) |
void | submit_gradient (const gradient_type &grad_in) |
template<typename Functor > | |
void | apply_quad_point_operations (const Functor &func) |
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 | |
types::global_dof_index * | local_to_global |
unsigned int | n_cells |
unsigned int | padding_length |
const unsigned int | mf_object_id |
const unsigned int | constraint_mask |
const bool | use_coloring |
Number * | inv_jac |
Number * | JxW |
Number * | values |
Number * | gradients [dim] |
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 90 of file cuda_fe_evaluation.h.
using CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::value_type = Number |
An alias for scalar quantities.
Definition at line 96 of file cuda_fe_evaluation.h.
using CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::gradient_type = Tensor<1, dim, Number> |
An alias for vectorial quantities.
Definition at line 101 of file cuda_fe_evaluation.h.
using CUDAWrappers::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 106 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 339 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 366 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 391 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 417 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 456 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 evaluate(true,...)
.
Definition at line 496 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 |
Same as above, except that the quadrature point is computed from thread id.
Definition at line 514 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_dof_value | ( | const unsigned int | dof | ) | const |
Return the value of a finite element function at degree of freedom dof
after a call to integrate() or before a call to evaluate().
Definition at line 533 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_dof_value |
Same as above, except that the local dof index is computed from the thread id.
Definition at line 551 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 566 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 | ) |
Same as above, except that the quadrature point is computed from the thread id.
Definition at line 580 of file cuda_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::submit_dof_value | ( | const value_type & | val_in, |
const unsigned int | dof | ||
) |
Write a value to the field containing the values for the degree of freedom with index dof
after a call to integrate() or before calling evaluate(). Access through the same fields as through get_dof_value().
Definition at line 595 of file cuda_fe_evaluation.h.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::submit_dof_value | ( | const value_type & | val_in | ) |
Same as above, except that the local dof index is computed from the thread id.
Definition at line 609 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 628 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 |
Same as above, except that the quadrature point is computed from the thread id.
Definition at line 660 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 690 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 | ) |
Same as above, except that the quadrature point is computed from the thread id.
Definition at line 713 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 738 of file cuda_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 755 of file cuda_fe_evaluation.h.
|
staticconstexpr |
Dimension.
Definition at line 111 of file cuda_fe_evaluation.h.
|
staticconstexpr |
Number of components.
Definition at line 116 of file cuda_fe_evaluation.h.
|
staticconstexpr |
Number of quadrature points per cell.
Definition at line 121 of file cuda_fe_evaluation.h.
|
staticconstexpr |
Number of tensor degrees of freedoms per cell.
Definition at line 127 of file cuda_fe_evaluation.h.
|
private |
Definition at line 314 of file cuda_fe_evaluation.h.
|
private |
Definition at line 315 of file cuda_fe_evaluation.h.
|
private |
Definition at line 316 of file cuda_fe_evaluation.h.
|
private |
Definition at line 317 of file cuda_fe_evaluation.h.
|
private |
Definition at line 319 of file cuda_fe_evaluation.h.
|
private |
Definition at line 321 of file cuda_fe_evaluation.h.
|
private |
Definition at line 323 of file cuda_fe_evaluation.h.
|
private |
Definition at line 324 of file cuda_fe_evaluation.h.
|
private |
Definition at line 327 of file cuda_fe_evaluation.h.
|
private |
Definition at line 328 of file cuda_fe_evaluation.h.