Reference documentation for deal.II version GIT relicensing-216-gcda843ca70 2024-03-28 12:20:02+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Attributes | List of all members
internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 > Struct Template Reference

#include <deal.II/matrix_free/tensor_product_kernels.h>

Public Member Functions

 EvaluatorTensorProduct ()
 
 EvaluatorTensorProduct (const AlignedVector< Number2 > &shape_values, const AlignedVector< Number2 > &shape_gradients, const AlignedVector< Number2 > &shape_hessians, const unsigned int=0, const unsigned int=0)
 
 EvaluatorTensorProduct (const Number2 *shape_values, const Number2 *shape_gradients, const Number2 *shape_hessians, const unsigned int dummy1=0, const unsigned int dummy2=0)
 
template<int direction, bool contract_over_rows, bool add, int stride = 1>
void values (const Number in[], Number out[]) const
 
template<int direction, bool contract_over_rows, bool add, int stride = 1>
void gradients (const Number in[], Number out[]) const
 
template<int direction, bool contract_over_rows, bool add>
void hessians (const Number in[], Number out[]) const
 
template<int direction, bool contract_over_rows, bool add>
void values_one_line (const Number in[], Number out[]) const
 
template<int direction, bool contract_over_rows, bool add>
void gradients_one_line (const Number in[], Number out[]) const
 
template<int direction, bool contract_over_rows, bool add>
void hessians_one_line (const Number in[], Number out[]) const
 

Static Public Member Functions

template<int direction, bool contract_over_rows, bool add, bool one_line = false, EvaluatorQuantity quantity = EvaluatorQuantity::value, int stride = 1>
static void apply (const Number2 *DEAL_II_RESTRICT shape_data, const Number *in, Number *out)
 

Static Public Attributes

static constexpr unsigned int n_rows_of_product
 
static constexpr unsigned int n_columns_of_product
 

Private Attributes

const Number2 * shape_values
 
const Number2 * shape_gradients
 
const Number2 * shape_hessians
 

Detailed Description

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number, typename Number2 = Number>
struct internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >

Generic evaluator framework that valuates the given shape data in general dimensions using the tensor product form. Depending on the particular layout in the matrix entries, this corresponds to a usual matrix-matrix product or a matrix-matrix product including some symmetries. The actual work is implemented by functions of type apply_matrix_vector_product working on a single dimension, controlled by suitable strides, using the kernel specified via variant.

Template Parameters
dimSpace dimension in which this class is applied
n_rowsNumber of rows in the transformation matrix, which corresponds to the number of 1d shape functions in the usual tensor contraction setting
n_columnsNumber of columns in the transformation matrix, which corresponds to the number of 1d shape functions in the usual tensor contraction setting
NumberAbstract number type for input and output arrays
Number2Abstract number type for coefficient arrays (defaults to same type as the input/output arrays); must implement operator* with Number and produce Number as an output to be a valid type

Definition at line 1137 of file tensor_product_kernels.h.

Constructor & Destructor Documentation

◆ EvaluatorTensorProduct() [1/3]

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::EvaluatorTensorProduct ( )
inline

Empty constructor. Does nothing. Be careful when using 'values' and related methods because they need to be filled with the other pointer

Definition at line 1148 of file tensor_product_kernels.h.

◆ EvaluatorTensorProduct() [2/3]

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::EvaluatorTensorProduct ( const AlignedVector< Number2 > &  shape_values,
const AlignedVector< Number2 > &  shape_gradients,
const AlignedVector< Number2 > &  shape_hessians,
const unsigned int  = 0,
const unsigned int  = 0 
)
inline

Constructor, taking the data from ShapeInfo

Definition at line 1157 of file tensor_product_kernels.h.

◆ EvaluatorTensorProduct() [3/3]

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::EvaluatorTensorProduct ( const Number2 *  shape_values,
const Number2 *  shape_gradients,
const Number2 *  shape_hessians,
const unsigned int  dummy1 = 0,
const unsigned int  dummy2 = 0 
)
inline

Constructor, taking the data from ShapeInfo via raw pointers

Definition at line 1197 of file tensor_product_kernels.h.

Member Function Documentation

◆ values()

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
template<int direction, bool contract_over_rows, bool add, int stride = 1>
void internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::values ( const Number  in[],
Number  out[] 
) const
inline

This interpolates values with sum factorization, according to the following parameters:

Template Parameters
directionThe direction along which the one-dimensional operations should be performed, using 0 as the fastest running direction x, 1 for y, and so on.
contract_over_rowsDescribes whether the interpolation is over the rows or the columns in the underlying interpolation matrix. With the chosen convention in the data fields, contract_over_rows==true means interpolation from DoF values to quadrature points, whereas using the argument false will perform summation over quadrature points to produce integrals for each test function.
addSpecify whether to add into the output array or overwrite the previous content.
strideThis parameter can specify an additional stride in the array associated with quadrature points (in for dof_to_quad==false, otherwise out) from consecutive points in the x-direction. This is used to place results from different interpolation steps next to each other in memory.
Parameters
inInput array for the operation, needs to be backed up by a sufficiently large memory region.
outArray holding the result of the sum factorization operation.

Definition at line 1237 of file tensor_product_kernels.h.

◆ gradients()

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
template<int direction, bool contract_over_rows, bool add, int stride = 1>
void internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::gradients ( const Number  in[],
Number  out[] 
) const
inline

This interpolates gradients with sum factorization, based on the second argument given to the constructor of this class. For the documentation of the template and function parameters, see the values function.

Definition at line 1251 of file tensor_product_kernels.h.

◆ hessians()

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
template<int direction, bool contract_over_rows, bool add>
void internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::hessians ( const Number  in[],
Number  out[] 
) const
inline

This interpolates hessians with sum factorization, based on the third argument given to the constructor of this class. For the documentation of the template and function parameters, see the values function.

Definition at line 1267 of file tensor_product_kernels.h.

◆ values_one_line()

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
template<int direction, bool contract_over_rows, bool add>
void internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::values_one_line ( const Number  in[],
Number  out[] 
) const
inline

A variant of interpolation with sum factorization that only applies the operation to a single 1d line, leaving all other entries untouched, rather than expanding the loop over all other directions of a tensor product mesh. For the documentation of the template and function parameters, see the other values function.

Definition at line 1287 of file tensor_product_kernels.h.

◆ gradients_one_line()

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
template<int direction, bool contract_over_rows, bool add>
void internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::gradients_one_line ( const Number  in[],
Number  out[] 
) const
inline

A variant of interpolation with sum factorization that only applies the gradient operation to a single 1d line, leaving all other entries untouched, rather than expanding the loop over all other directions of a tensor product mesh. For the documentation of the template and function parameters, see the values function.

Definition at line 1303 of file tensor_product_kernels.h.

◆ hessians_one_line()

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
template<int direction, bool contract_over_rows, bool add>
void internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::hessians_one_line ( const Number  in[],
Number  out[] 
) const
inline

A variant of interpolation with sum factorization that only applies the Hessian operation to a single 1d line, leaving all other entries untouched, rather than expanding the loop over all other directions of a tensor product mesh. For the documentation of the template and function parameters, see the values function.

Definition at line 1322 of file tensor_product_kernels.h.

◆ apply()

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 >
template<int direction, bool contract_over_rows, bool add, bool one_line, EvaluatorQuantity quantity, int stride>
void internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::apply ( const Number2 *DEAL_II_RESTRICT  shape_data,
const Number *  in,
Number *  out 
)
inlinestatic

This function applies the tensor product kernel with sum factorization, corresponding to a matrix-vector multiplication of 1d stripes, along the given direction of the tensor data in the input array. This function allows the in and out arrays to alias for the case n_rows == n_columns, i.e., it is safe to perform the contraction in place where in and out point to the same address. For the case n_rows != n_columns, the output is in general not correct.

Template Parameters
directionDirection that is evaluated
contract_over_rowsIf true, the tensor contraction sums over the rows in the given shape_data array, otherwise it sums over the columns
addIf true, the result is added to the output vector, else the computed values overwrite the content in the output
one_lineIf true, the kernel is only applied along a single 1d stripe within a dim-dimensional tensor, not the full n_rows^dim points as in the false case.
quantitySpecify whether values, gradients or Hessians should be interpolated, allowing specialized algorithms for some class template parameters of variant to find the right path.
strideThis parameter enables to place the result of the tensor product evaluation in the output array (if contract_over_rows == true) or input array (if contract_over_rows == false) with additional strides between adjacent points in x direction, which is used to group all components of a gradient adjacent in memory. If the stride is one, the data will form a contiguous range in memory.
Parameters
shape_dataTransformation matrix with n_rows rows and n_columns columns, stored in row-major format.
inPointer to the start of the input data vector.
outPointer to the start of the output data vector.

Definition at line 1402 of file tensor_product_kernels.h.

Member Data Documentation

◆ n_rows_of_product

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
constexpr unsigned int internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::n_rows_of_product
staticconstexpr
Initial value:
=
Utilities::pow(n_rows, dim)
constexpr T pow(const T base, const int iexp)
Definition utilities.h:963

Definition at line 1139 of file tensor_product_kernels.h.

◆ n_columns_of_product

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
constexpr unsigned int internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::n_columns_of_product
staticconstexpr
Initial value:
=
Utilities::pow(n_columns, dim)

Definition at line 1141 of file tensor_product_kernels.h.

◆ shape_values

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
const Number2* internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::shape_values
private

Definition at line 1382 of file tensor_product_kernels.h.

◆ shape_gradients

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
const Number2* internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::shape_gradients
private

Definition at line 1383 of file tensor_product_kernels.h.

◆ shape_hessians

template<EvaluatorVariant variant, int dim, int n_rows, int n_columns, typename Number , typename Number2 = Number>
const Number2* internal::EvaluatorTensorProduct< variant, dim, n_rows, n_columns, Number, Number2 >::shape_hessians
private

Definition at line 1384 of file tensor_product_kernels.h.


The documentation for this struct was generated from the following file: