Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Classes | Typedefs | Functions
Differentiation::AD::internal Namespace Reference

Classes

struct  ADNumberInfoFromEnum
 
struct  ExtractData
 
struct  Extractor
 
struct  Extractor< dim, FEValuesExtractors::Scalar >
 
struct  Extractor< dim, FEValuesExtractors::SymmetricTensor< 2 > >
 
struct  Extractor< dim, FEValuesExtractors::Tensor< 1 > >
 
struct  Extractor< dim, FEValuesExtractors::Tensor< 2 > >
 
struct  Extractor< dim, FEValuesExtractors::Vector >
 
struct  HasRequiredADInfo
 
struct  HessianType
 
struct  HessianType< FEValuesExtractors::Scalar, FEValuesExtractors::SymmetricTensor< 2 > >
 
struct  HessianType< FEValuesExtractors::SymmetricTensor< 2 >, FEValuesExtractors::Scalar >
 
struct  HessianType< FEValuesExtractors::SymmetricTensor< 2 >, FEValuesExtractors::SymmetricTensor< 2 > >
 
struct  Marking
 
struct  NumberType
 
struct  RemoveComplexWrapper
 
struct  SacadoNumberInfo
 
struct  ScalarFieldGradient
 
struct  ScalarFieldHessian
 

Typedefs

template<int dim, typename NumberType , typename ExtractorType_Field >
using VectorFieldValue = ScalarFieldGradient< dim, NumberType, ExtractorType_Field >
 
template<int dim, typename NumberType , typename ExtractorType_Field , typename ExtractorType_Derivative >
using VectorFieldJacobian = ScalarFieldHessian< dim, NumberType, ExtractorType_Field, ExtractorType_Derivative >
 

Functions

template<int dim, typename IndexType = unsigned int, typename ExtractorType >
std::vector< IndexType > extract_field_component_indices (const ExtractorType &extractor, const bool ignore_symmetries=true)
 
template<int dim, typename IndexType = unsigned int>
std::vector< IndexType > extract_field_component_indices (const FEValuesExtractors::SymmetricTensor< 2 > &extractor_symm_tensor, const bool ignore_symmetries=true)
 
template<typename TensorType , typename NumberType >
void set_tensor_entry (TensorType &t, const unsigned int &unrolled_index, const NumberType &value)
 
template<int dim, typename NumberType >
void set_tensor_entry (Tensor< 0, dim, NumberType > &t, const unsigned int &unrolled_index, const NumberType &value)
 
template<typename NumberType >
void set_tensor_entry (NumberType &t, const unsigned int &unrolled_index, const NumberType &value)
 
template<int dim, typename NumberType >
void set_tensor_entry (SymmetricTensor< 4, dim, NumberType > &t, const unsigned int &unrolled_index_row, const unsigned int &unrolled_index_col, const NumberType &value)
 
template<int rank, int dim, typename NumberType , template< int, int, typename > class TensorType>
NumberType get_tensor_entry (const TensorType< rank, dim, NumberType > &t, const unsigned int &unrolled_index)
 
template<int dim, typename NumberType , template< int, int, typename > class TensorType>
NumberType get_tensor_entry (const TensorType< 0, dim, NumberType > &t, const unsigned int &unrolled_index)
 
template<typename NumberType >
const NumberTypeget_tensor_entry (const NumberType &t, const unsigned int &unrolled_index)
 
template<int rank, int dim, typename NumberType , template< int, int, typename > class TensorType>
NumberTypeget_tensor_entry (TensorType< rank, dim, NumberType > &t, const unsigned int &unrolled_index)
 
template<int dim, typename NumberType , template< int, int, typename > class TensorType>
NumberTypeget_tensor_entry (TensorType< 0, dim, NumberType > &t, const unsigned int &index)
 
template<typename NumberType >
NumberTypeget_tensor_entry (NumberType &t, const unsigned int &index)
 
template<typename ADNumberType >
std::enable_if<!(ADNumberTraits< ADNumberType >::type_code==NumberTypes::sacado_rad||ADNumberTraits< ADNumberType >::type_code==NumberTypes::sacado_rad_dfad)>::type reverse_mode_dependent_variable_activation (ADNumberType &)
 
template<typename ADNumberType >
std::enable_if< ADNumberTraits< ADNumberType >::type_code==NumberTypes::sacado_rad||ADNumberTraits< ADNumberType >::type_code==NumberTypes::sacado_rad_dfad >::type reverse_mode_dependent_variable_activation (ADNumberType &dependent_variable)
 
template<typename ADNumberType >
std::enable_if<!(ADNumberTraits< ADNumberType >::type_code==NumberTypes::adolc_tapeless)>::type configure_tapeless_mode (const unsigned int)
 
template<typename ADNumberType >
std::enable_if< ADNumberTraits< ADNumberType >::type_code==NumberTypes::adolc_tapeless >::type configure_tapeless_mode (const unsigned int n_directional_derivatives)
 

Detailed Description

This namespace defines the classes that help provide a unified interface to all auto-differentiation numbers.

Typedef Documentation

◆ VectorFieldValue

template<int dim, typename NumberType , typename ExtractorType_Field >
using Differentiation::AD::internal::VectorFieldValue = typedef ScalarFieldGradient<dim, NumberType, ExtractorType_Field>

A helper struct that defines the return type of value computations of vector fields the ExtractorType_Field template parameter.

Definition at line 2336 of file ad_helpers.h.

◆ VectorFieldJacobian

template<int dim, typename NumberType , typename ExtractorType_Field , typename ExtractorType_Derivative >
using Differentiation::AD::internal::VectorFieldJacobian = typedef ScalarFieldHessian<dim, NumberType, ExtractorType_Field, ExtractorType_Derivative>

A helper struct that defines the final return type of Jacobian (first derivative) calculations of vector fields with respect to fields as defined by the two extractor-type template parameters. The first, ExtractorType_Field, defines the field from which the initial field values are computed while the second, ExtractorType_Derivative, defines the field that the derivatives are taken with respect to.

Definition at line 2355 of file ad_helpers.h.

Function Documentation

◆ extract_field_component_indices() [1/2]

template<int dim, typename IndexType = unsigned int, typename ExtractorType >
std::vector<IndexType> Differentiation::AD::internal::extract_field_component_indices ( const ExtractorType &  extractor,
const bool  ignore_symmetries = true 
)

Return a global view of the field component indices that correspond to the input extractor. For this general function the ignore_symmetries flag has no effect.

Definition at line 2367 of file ad_helpers.h.

◆ extract_field_component_indices() [2/2]

template<int dim, typename IndexType = unsigned int>
std::vector<IndexType> Differentiation::AD::internal::extract_field_component_indices ( const FEValuesExtractors::SymmetricTensor< 2 > &  extractor_symm_tensor,
const bool  ignore_symmetries = true 
)

Return a global view of the field component indices that correspond to the input FEValuesExtractors::SymmetricTensor extractor_symm_tensor. If the ignore_symmetries is set true, then all component of the tensor are considered to be independent. If set to false, then the set of returned indices will contain duplicate entries for components that are symmetric.

Definition at line 2391 of file ad_helpers.h.

◆ set_tensor_entry() [1/4]

template<typename TensorType , typename NumberType >
void Differentiation::AD::internal::set_tensor_entry ( TensorType &  t,
const unsigned int &  unrolled_index,
const NumberType value 
)
inline

Set the unrolled component given by index in the generic tensor t to the given value.

Definition at line 2441 of file ad_helpers.h.

◆ set_tensor_entry() [2/4]

template<int dim, typename NumberType >
void Differentiation::AD::internal::set_tensor_entry ( Tensor< 0, dim, NumberType > &  t,
const unsigned int &  unrolled_index,
const NumberType value 
)
inline

Set the unrolled component given by index in the rank-0 tensor t to the given value.

Definition at line 2457 of file ad_helpers.h.

◆ set_tensor_entry() [3/4]

template<typename NumberType >
void Differentiation::AD::internal::set_tensor_entry ( NumberType t,
const unsigned int &  unrolled_index,
const NumberType value 
)
inline

Set the value of t to the given value. This function exists to provide compatibility with similar functions that exist for use with the tensor classes.

Definition at line 2474 of file ad_helpers.h.

◆ set_tensor_entry() [4/4]

template<int dim, typename NumberType >
void Differentiation::AD::internal::set_tensor_entry ( SymmetricTensor< 4, dim, NumberType > &  t,
const unsigned int &  unrolled_index_row,
const unsigned int &  unrolled_index_col,
const NumberType value 
)
inline

Set the unrolled component given by the index_row and the index_col in the fourth-order symmetric tensor t to the given value.

Definition at line 2490 of file ad_helpers.h.

◆ get_tensor_entry() [1/6]

template<int rank, int dim, typename NumberType , template< int, int, typename > class TensorType>
NumberType Differentiation::AD::internal::get_tensor_entry ( const TensorType< rank, dim, NumberType > &  t,
const unsigned int &  unrolled_index 
)
inline

Return the value of the index'th unrolled component of the generic tensor t.

Definition at line 2524 of file ad_helpers.h.

◆ get_tensor_entry() [2/6]

template<int dim, typename NumberType , template< int, int, typename > class TensorType>
NumberType Differentiation::AD::internal::get_tensor_entry ( const TensorType< 0, dim, NumberType > &  t,
const unsigned int &  unrolled_index 
)
inline

Return the value of the index'th unrolled component of the rank-0 tensor t.

Definition at line 2543 of file ad_helpers.h.

◆ get_tensor_entry() [3/6]

template<typename NumberType >
const NumberType& Differentiation::AD::internal::get_tensor_entry ( const NumberType t,
const unsigned int &  unrolled_index 
)
inline

Return the value of t. This function exists to provide compatibility with similar functions that exist for use with the tensor classes.

Definition at line 2559 of file ad_helpers.h.

◆ get_tensor_entry() [4/6]

template<int rank, int dim, typename NumberType , template< int, int, typename > class TensorType>
NumberType& Differentiation::AD::internal::get_tensor_entry ( TensorType< rank, dim, NumberType > &  t,
const unsigned int &  unrolled_index 
)
inline

Return a reference to the entry stored in the index'th unrolled component of the generic tensor t.

Definition at line 2576 of file ad_helpers.h.

◆ get_tensor_entry() [5/6]

template<int dim, typename NumberType , template< int, int, typename > class TensorType>
NumberType& Differentiation::AD::internal::get_tensor_entry ( TensorType< 0, dim, NumberType > &  t,
const unsigned int &  index 
)

Return a reference to the entry stored in the index'th unrolled component of the rank-0 tensor t.

Definition at line 2594 of file ad_helpers.h.

◆ get_tensor_entry() [6/6]

template<typename NumberType >
NumberType& Differentiation::AD::internal::get_tensor_entry ( NumberType t,
const unsigned int &  index 
)
inline

Return a reference to t. This function exists to provide compatibility with similar functions that exist for use with the tensor classes.

Definition at line 2610 of file ad_helpers.h.

◆ reverse_mode_dependent_variable_activation() [1/2]

template<typename ADNumberType >
std::enable_if<!(ADNumberTraits<ADNumberType>::type_code == NumberTypes::sacado_rad || ADNumberTraits<ADNumberType>::type_code == NumberTypes::sacado_rad_dfad)>::type Differentiation::AD::internal::reverse_mode_dependent_variable_activation ( ADNumberType &  )

A dummy function to define the active dependent variable when using reverse-mode AD.

Definition at line 1583 of file ad_drivers.cc.

◆ reverse_mode_dependent_variable_activation() [2/2]

template<typename ADNumberType >
std::enable_if<ADNumberTraits<ADNumberType>::type_code == NumberTypes::sacado_rad || ADNumberTraits<ADNumberType>::type_code == NumberTypes::sacado_rad_dfad>::type Differentiation::AD::internal::reverse_mode_dependent_variable_activation ( ADNumberType &  dependent_variable)

Define the active dependent variable when using reverse-mode AD.

If there are multiple dependent variables then it is necessary to inform the independent variables, from which the adjoints are computed, which dependent variable they are computing the gradients with respect to. This function broadcasts this information.

Definition at line 1602 of file ad_drivers.cc.

◆ configure_tapeless_mode() [1/2]

template<typename ADNumberType >
std::enable_if<!(ADNumberTraits<ADNumberType>::type_code == NumberTypes::adolc_tapeless)>::type Differentiation::AD::internal::configure_tapeless_mode ( const unsigned int  )

A dummy function to enable vector mode for tapeless auto-differentiable numbers.

Definition at line 1626 of file ad_drivers.cc.

◆ configure_tapeless_mode() [2/2]

template<typename ADNumberType >
std::enable_if<ADNumberTraits<ADNumberType>::type_code == NumberTypes::adolc_tapeless>::type Differentiation::AD::internal::configure_tapeless_mode ( const unsigned int  n_directional_derivatives)

Enable vector mode for ADOL-C tapeless numbers.

This function checks to see if its legal to increase the maximum number of directional derivatives to be considered during calculations. If not then it throws an error.

Definition at line 1643 of file ad_drivers.cc.