Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Public Types | Public Member Functions | Static Public Attributes | Private Attributes | List of all members
CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number > Class Template Reference

#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_indexlocal_to_global
 
unsigned int n_cells
 
unsigned int padding_length
 
const unsigned int constraint_mask
 
const bool use_coloring
 
Number * inv_jac
 
Number * JxW
 
Number * values
 
Number * gradients [dim]
 

Detailed Description

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
class CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >

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:

Template Parameters
dimDimension in which this class is to be used
fe_degreeDegree of the tensor prodict finite element with fe_degree+1 degrees of freedom per coordinate direction
n_q_points_1dNumber of points in the quadrature formular in 1D, defaults to fe_degree+1
n_componentsNumber 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
NumberNumber format, double or float. Defaults to double.
Author
Karl Ljungkvist, Bruno Turcksin, 2016

Definition at line 92 of file cuda_fe_evaluation.h.

Member Typedef Documentation

◆ value_type

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
using CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::value_type = Number

An alias for scalar quantities.

Definition at line 98 of file cuda_fe_evaluation.h.

◆ gradient_type

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
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 103 of file cuda_fe_evaluation.h.

◆ data_type

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
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 108 of file cuda_fe_evaluation.h.

Constructor & Destructor Documentation

◆ FEEvaluation()

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 341 of file cuda_fe_evaluation.h.

Member Function Documentation

◆ read_dof_values()

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 367 of file cuda_fe_evaluation.h.

◆ distribute_local_to_global()

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 392 of file cuda_fe_evaluation.h.

◆ evaluate()

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 418 of file cuda_fe_evaluation.h.

◆ integrate()

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 457 of file cuda_fe_evaluation.h.

◆ get_value() [1/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 497 of file cuda_fe_evaluation.h.

◆ get_value() [2/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 516 of file cuda_fe_evaluation.h.

◆ get_dof_value() [1/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 535 of file cuda_fe_evaluation.h.

◆ get_dof_value() [2/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 553 of file cuda_fe_evaluation.h.

◆ submit_value() [1/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 568 of file cuda_fe_evaluation.h.

◆ submit_value() [2/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 582 of file cuda_fe_evaluation.h.

◆ submit_dof_value() [1/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 597 of file cuda_fe_evaluation.h.

◆ submit_dof_value() [2/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 611 of file cuda_fe_evaluation.h.

◆ get_gradient() [1/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 630 of file cuda_fe_evaluation.h.

◆ get_gradient() [2/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 662 of file cuda_fe_evaluation.h.

◆ submit_gradient() [1/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 692 of file cuda_fe_evaluation.h.

◆ submit_gradient() [2/2]

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
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 715 of file cuda_fe_evaluation.h.

◆ apply_quad_point_operations()

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
template<typename Functor >
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

__device__ void operator()(
const unsigned int q_point) const;

Definition at line 740 of file cuda_fe_evaluation.h.

◆ apply_for_each_quad_point()

template<int dim, int fe_degree, int n_q_points_1d, int n_components_, typename Number >
template<typename Functor >
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 757 of file cuda_fe_evaluation.h.

Member Data Documentation

◆ dimension

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
constexpr unsigned int CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::dimension = dim
staticconstexpr

Dimension.

Definition at line 113 of file cuda_fe_evaluation.h.

◆ n_components

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
constexpr unsigned int CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::n_components = n_components_
staticconstexpr

Number of components.

Definition at line 118 of file cuda_fe_evaluation.h.

◆ n_q_points

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
constexpr unsigned int CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::n_q_points
staticconstexpr
Initial value:
=
Utilities::pow(n_q_points_1d, dim)

Number of quadrature points per cell.

Definition at line 123 of file cuda_fe_evaluation.h.

◆ tensor_dofs_per_cell

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
constexpr unsigned int CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::tensor_dofs_per_cell
staticconstexpr
Initial value:
=
Utilities::pow(fe_degree + 1, dim)

Number of tensor degrees of freedoms per cell.

Definition at line 129 of file cuda_fe_evaluation.h.

◆ local_to_global

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
types::global_dof_index* CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::local_to_global
private

Definition at line 316 of file cuda_fe_evaluation.h.

◆ n_cells

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
unsigned int CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::n_cells
private

Definition at line 317 of file cuda_fe_evaluation.h.

◆ padding_length

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
unsigned int CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::padding_length
private

Definition at line 318 of file cuda_fe_evaluation.h.

◆ constraint_mask

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
const unsigned int CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::constraint_mask
private

Definition at line 320 of file cuda_fe_evaluation.h.

◆ use_coloring

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
const bool CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::use_coloring
private

Definition at line 322 of file cuda_fe_evaluation.h.

◆ inv_jac

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
Number* CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::inv_jac
private

Definition at line 324 of file cuda_fe_evaluation.h.

◆ JxW

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
Number* CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::JxW
private

Definition at line 325 of file cuda_fe_evaluation.h.

◆ values

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
Number* CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::values
private

Definition at line 328 of file cuda_fe_evaluation.h.

◆ gradients

template<int dim, int fe_degree, int n_q_points_1d = fe_degree + 1, int n_components_ = 1, typename Number = double>
Number* CUDAWrappers::FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number >::gradients[dim]
private

Definition at line 329 of file cuda_fe_evaluation.h.


The documentation for this class was generated from the following file:
Utilities::pow
constexpr T pow(const T base, const int iexp)
Definition: utilities.h:476
CUDAWrappers::FEEvaluation
Definition: cuda_fe_evaluation.h:92