Reference documentation for deal.II version 9.5.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\}}\)
Loading...
Searching...
No Matches
List of all members
Differentiation::AD::TapedDrivers< ADNumberType, ScalarType, T > Struct Template Reference

#include <deal.II/differentiation/ad/ad_drivers.h>

Inheritance diagram for Differentiation::AD::TapedDrivers< ADNumberType, ScalarType, T >:
[legend]

Public Member Functions

Taping
bool is_recording () const
 
Types< ADNumberType >::tape_index active_tape_index () const
 
bool is_registered_tape (const typename Types< ADNumberType >::tape_index tape_index) const
 
bool keep_independent_values () const
 
void set_tape_buffer_sizes (const typename Types< ADNumberType >::tape_buffer_sizes obufsize=64 *1024 *1024, const typename Types< ADNumberType >::tape_buffer_sizes lbufsize=64 *1024 *1024, const typename Types< ADNumberType >::tape_buffer_sizes vbufsize=64 *1024 *1024, const typename Types< ADNumberType >::tape_buffer_sizes tbufsize=64 *1024 *1024)
 
void start_taping (const typename Types< ADNumberType >::tape_index tape_index, const bool keep_independent_values)
 
void stop_taping (const typename Types< ADNumberType >::tape_index active_tape_index, const bool write_tapes_to_file)
 
std::vector< typename Types< ADNumberType >::tape_index > get_registered_tape_indices () const
 
void activate_tape (const typename Types< ADNumberType >::tape_index tape_index)
 
bool requires_retaping (const typename Types< ADNumberType >::tape_index tape_index) const
 
bool last_action_requires_retaping () const
 
void remove_tape (const typename Types< ADNumberType >::tape_index tape_index)
 
void reset (const bool clear_registered_tapes)
 
void print (std::ostream &stream) const
 
void print_tape_stats (const typename Types< ADNumberType >::tape_index tape_index, std::ostream &stream) const
 
Drivers for scalar functions (one dependent variable)
ScalarType value (const typename Types< ADNumberType >::tape_index active_tape_index, const std::vector< ScalarType > &independent_variables) const
 
void gradient (const typename Types< ADNumberType >::tape_index active_tape_index, const std::vector< ScalarType > &independent_variables, Vector< ScalarType > &gradient) const
 
void hessian (const typename Types< ADNumberType >::tape_index active_tape_index, const std::vector< ScalarType > &independent_variables, FullMatrix< ScalarType > &hessian) const
 
Drivers for vector functions (multiple dependent variables)
void values (const typename Types< ADNumberType >::tape_index active_tape_index, const unsigned int n_dependent_variables, const std::vector< ScalarType > &independent_variables, Vector< ScalarType > &values) const
 
void jacobian (const typename Types< ADNumberType >::tape_index active_tape_index, const unsigned int n_dependent_variables, const std::vector< ScalarType > &independent_variables, FullMatrix< ScalarType > &jacobian) const
 

Detailed Description

template<typename ADNumberType, typename ScalarType, typename T = void>
struct Differentiation::AD::TapedDrivers< ADNumberType, ScalarType, T >

A prototype driver class for taped auto-differentiable numbers.

It is intended that this class be specialized for the valid combinations of auto-differentiable numbers and output scalar number types.

Template Parameters
ADNumberTypeA type corresponding to a supported auto-differentiable number.
ScalarTypeA real or complex floating point number type that is the scalar value type used for input to, and output from, operations performed with auto-differentiable numbers.
TAn arbitrary type resulting from the application of the SFINAE idiom to selectively specialize this class.

Definition at line 147 of file ad_drivers.h.

Member Function Documentation

◆ is_recording()

template<typename ADNumberType >
bool Differentiation::AD::TapedDrivers< ADNumberType >::is_recording

Return whether or not this class is tracking calculations performed with its marked independent variables.

Definition at line 54 of file ad_drivers.cc.

◆ active_tape_index()

template<typename ADNumberType >
Types< ADNumberType >::tape_index Differentiation::AD::TapedDrivers< ADNumberType >::active_tape_index

Return the tape number which is currently activated for recording or reading.

Definition at line 63 of file ad_drivers.cc.

◆ is_registered_tape()

template<typename ADNumberType >
bool Differentiation::AD::TapedDrivers< ADNumberType >::is_registered_tape ( const typename Types< ADNumberType >::tape_index  tape_index) const

Return whether or not a tape number has already been used or registered.

Definition at line 72 of file ad_drivers.cc.

◆ keep_independent_values()

template<typename ADNumberType >
bool Differentiation::AD::TapedDrivers< ADNumberType >::keep_independent_values

Return whether or not the numerical values of all independent variables are recorded in the tape buffer.

Definition at line 82 of file ad_drivers.cc.

◆ set_tape_buffer_sizes()

template<typename ADNumberType >
void Differentiation::AD::TapedDrivers< ADNumberType >::set_tape_buffer_sizes ( const typename Types< ADNumberType >::tape_buffer_sizes  obufsize = 64 * 1024 * 1024,
const typename Types< ADNumberType >::tape_buffer_sizes  lbufsize = 64 * 1024 * 1024,
const typename Types< ADNumberType >::tape_buffer_sizes  vbufsize = 64 * 1024 * 1024,
const typename Types< ADNumberType >::tape_buffer_sizes  tbufsize = 64 * 1024 * 1024 
)

Set the buffer sizes for the next active tape.

This function must be called before start_recording_operations() for it to have any influence on the memory allocated to the next recorded tape.

Note
This function only has an effect when using ADOL-C numbers. As stated by the ADOL-C manual, it may be desirable to create a file ".adolcrc" in the program run directory and set the buffer size therein. Alternatively, this function can be used to override the settings for any given tape, or can be used in the event that no ".adolcrc" file exists. The default value for each buffer is set at 64 megabytes, a heuristically chosen value thought to be appropriate for use within the context of finite element analysis when considering coupled problems with multiple vector-valued fields discretised by higher order shape functions, as well as complex constitutive laws.
Parameters
[in]obufsizeADOL-C operations buffer size
[in]lbufsizeADOL-C locations buffer size
[in]vbufsizeADOL-C value buffer size
[in]tbufsizeADOL-C Taylor buffer size

Definition at line 91 of file ad_drivers.cc.

◆ start_taping()

template<typename ADNumberType >
void Differentiation::AD::TapedDrivers< ADNumberType >::start_taping ( const typename Types< ADNumberType >::tape_index  tape_index,
const bool  keep_independent_values 
)

Enable the recording mode for a given tape.

Parameters
[in]tape_indexThe index of the tape to be written.
[in]keep_independent_valuesDetermines whether the numerical values of all independent variables are recorded in the tape buffer.

Definition at line 103 of file ad_drivers.cc.

◆ stop_taping()

template<typename ADNumberType >
void Differentiation::AD::TapedDrivers< ADNumberType >::stop_taping ( const typename Types< ADNumberType >::tape_index  active_tape_index,
const bool  write_tapes_to_file 
)

Disable the recording mode for a given tape.

Parameters
[in]active_tape_indexThe index of the (currently active) tape to be finalized and potentially written to file.
[in]write_tapes_to_fileA flag that specified whether the tape should be written to file or kept in memory.

Definition at line 113 of file ad_drivers.cc.

◆ get_registered_tape_indices()

template<typename ADNumberType >
std::vector< typename Types< ADNumberType >::tape_index > Differentiation::AD::TapedDrivers< ADNumberType >::get_registered_tape_indices

Return a list of registered tape indices.

Definition at line 123 of file ad_drivers.cc.

◆ activate_tape()

template<typename ADNumberType >
void Differentiation::AD::TapedDrivers< ADNumberType >::activate_tape ( const typename Types< ADNumberType >::tape_index  tape_index)

Select a tape to record to or read from.

This function activates a tape, but depending on whether read_mode is set, the tape is either taken as previously written to (and put into read-only mode), or cleared for (re-)taping.

Parameters
[in]tape_indexThe index of the tape to be written to/read from.
Note
The chosen tape index must be greater than Numbers<ADNumberType>::invalid_tape_index and less than Numbers<ADNumberType>::max_tape_index.

Definition at line 133 of file ad_drivers.cc.

◆ requires_retaping()

template<typename ADNumberType >
bool Differentiation::AD::TapedDrivers< ADNumberType >::requires_retaping ( const typename Types< ADNumberType >::tape_index  tape_index) const

Return a flag that, when true, indicates that the retaping of the dependent function for the chosen tape_index is necessary for a reliable computation to be performed. This may be necessary if a sign comparison within branched operations yields different results to those computed at the original tape evaluation point.

This issue, known as "branch switching", can be clarified by means of a trivial, contrived example:

ADNumberType func (ADNumberType x, ADNumberType y, ADNumberType z)
{
if (x < y)
return y;
else
return x*z;
}

During taping, the conditional statement may be either true or false, and the result (with its sensitivities) returned by this function. The AD library doesn't just record the parse tree of the operations applied on the branch chosen at the time to taping, but also checks that the condition continues to be satisfied. For some other evaluation of the tape (i.e. for some different inputs x and y), the other branch of the conditional check may be chosen. The result of following this code path has not been recorded on the tape, and therefore cannot be evaluated. In such a case, the underlying AD library will be able to tell you that it is necessary to re-record the tape at the new evaluation point in order to resolve the new code branch. This function can be used to find out whether this is so.

Note
The chosen tape index must be greater than Numbers<ADNumberType>::invalid_tape_index and less than Numbers<ADNumberType>::max_tape_index.

Definition at line 142 of file ad_drivers.cc.

◆ last_action_requires_retaping()

template<typename ADNumberType >
bool Differentiation::AD::TapedDrivers< ADNumberType >::last_action_requires_retaping

Return a flag that, when true, indicates that the retaping of the dependent function is necessary for a reliable computation to be performed on the currently active tape. This may be necessary if a sign comparison within branched operations yields different results to those computed at the original tape evaluation point.

This issue, known as "branch switching", can be clarified by means of a trivial, contrived example:

ADNumberType func (ADNumberType x, ADNumberType y, ADNumberType z)
{
if (x < y)
return y;
else
return x*z;
}

During taping, the conditional statement may be either true or false, and the result (with its sensitivities) returned by this function. The AD library doesn't just record the parse tree of the operations applied on the branch chosen at the time to taping, but also checks that the condition continues to be satisfied. For some other evaluation of the tape (i.e. for some different inputs x and y), the other branch of the conditional check may be chosen. The result of following this code path has not been recorded on the tape, and therefore cannot be evaluated. In such a case, the underlying AD library will be able to tell you that it is necessary to re-record the tape at the new evaluation point in order to resolve the new code branch. This function can be used to find out whether this is so.

Definition at line 152 of file ad_drivers.cc.

◆ remove_tape()

template<typename ADNumberType >
void Differentiation::AD::TapedDrivers< ADNumberType >::remove_tape ( const typename Types< ADNumberType >::tape_index  tape_index)

Completely erases the tape with the given tape_index.

Definition at line 162 of file ad_drivers.cc.

◆ reset()

template<typename ADNumberType >
void Differentiation::AD::TapedDrivers< ADNumberType >::reset ( const bool  clear_registered_tapes)

Reset the state of the class.

Note
This also resets the active tape number to an invalid number, and deactivates the recording mode for taped variables.

Definition at line 171 of file ad_drivers.cc.

◆ print()

template<typename ADNumberType >
void Differentiation::AD::TapedDrivers< ADNumberType >::print ( std::ostream &  stream) const

Print the status of all queryable data. Exactly what is printed and its format depends on the ad_type, as is determined by the ADNumberTypeCode template parameter.

Parameters
[in]streamThe output stream to which the values are to be written.

Definition at line 179 of file ad_drivers.cc.

◆ print_tape_stats()

template<typename ADNumberType >
void Differentiation::AD::TapedDrivers< ADNumberType >::print_tape_stats ( const typename Types< ADNumberType >::tape_index  tape_index,
std::ostream &  stream 
) const

Print the statistics regarding the usage of the tapes.

Parameters
[in]tape_indexThe index of the tape to get the statistics of.
[out]streamThe output stream to which the values are to be written.

Definition at line 187 of file ad_drivers.cc.

◆ value()

template<typename ADNumberType , typename ScalarType , typename T >
TapedDrivers< ADNumberType, float, std::enable_if_t< ADNumberTraits< ADNumberType >::type_code==NumberTypes::adolc_taped > >::scalar_type Differentiation::AD::TapedDrivers< ADNumberType >::value ( const typename Types< ADNumberType >::tape_index  active_tape_index,
const std::vector< ScalarType > &  independent_variables 
) const

Compute the value of the scalar field.

Parameters
[in]active_tape_indexThe index of the tape on which the dependent function is recorded.
[in]independent_variablesThe scalar values of the independent variables whose sensitivities were tracked.
Returns
The scalar value of the function.

Definition at line 197 of file ad_drivers.cc.

◆ gradient()

template<typename ADNumberType , typename ScalarType , typename T >
void Differentiation::AD::TapedDrivers< ADNumberType >::gradient ( const typename Types< ADNumberType >::tape_index  active_tape_index,
const std::vector< ScalarType > &  independent_variables,
Vector< ScalarType > &  gradient 
) const

Compute the gradient of the scalar field with respect to all independent variables.

Parameters
[in]active_tape_indexThe index of the tape on which the dependent function is recorded.
[in]independent_variablesThe scalar values of the independent variables whose sensitivities were tracked.
[out]gradientThe values of the dependent function's gradients. It is expected that this vector be of the correct size (with length n_independent_variables).

Definition at line 208 of file ad_drivers.cc.

◆ hessian()

template<typename ADNumberType , typename ScalarType , typename T >
void Differentiation::AD::TapedDrivers< ADNumberType >::hessian ( const typename Types< ADNumberType >::tape_index  active_tape_index,
const std::vector< ScalarType > &  independent_variables,
FullMatrix< ScalarType > &  hessian 
) const

Compute the Hessian of the scalar field with respect to all independent variables.

Parameters
[in]active_tape_indexThe index of the tape on which the dependent function is recorded.
[in]independent_variablesThe scalar values of the independent variables whose sensitivities were tracked.
[out]hessianThe values of the dependent function's Hessian. It is expected that this matrix be of the correct size (with dimensions n_independent_variables \(\times\)n_independent_variables).

Definition at line 219 of file ad_drivers.cc.

◆ values()

template<typename ADNumberType , typename ScalarType , typename T >
void Differentiation::AD::TapedDrivers< ADNumberType >::values ( const typename Types< ADNumberType >::tape_index  active_tape_index,
const unsigned int  n_dependent_variables,
const std::vector< ScalarType > &  independent_variables,
Vector< ScalarType > &  values 
) const

Compute the values of the vector field.

Parameters
[in]active_tape_indexThe index of the tape on which the dependent function is recorded.
[in]n_dependent_variablesThe number of dependent variables.
[in]independent_variablesThe scalar values of the independent variables whose sensitivities were tracked.
[out]valuesThe component values of the dependent functions. It is expected that this vector be of the correct size (with length n_dependent_variables).

Definition at line 230 of file ad_drivers.cc.

◆ jacobian()

template<typename ADNumberType , typename ScalarType , typename T >
void Differentiation::AD::TapedDrivers< ADNumberType >::jacobian ( const typename Types< ADNumberType >::tape_index  active_tape_index,
const unsigned int  n_dependent_variables,
const std::vector< ScalarType > &  independent_variables,
FullMatrix< ScalarType > &  jacobian 
) const

Compute the Jacobian of the vector field.

The Jacobian of a vector field is in essence the gradient of each dependent variable with respect to all independent variables. This operation is therefore analogous to the gradient() operation performed on a collection of scalar valued fields.

Parameters
[in]active_tape_indexThe index of the tape on which the dependent function is recorded.
[in]n_dependent_variablesThe number of dependent variables.
[in]independent_variablesThe scalar values of the independent variables whose sensitivities were tracked.
[out]jacobianThe component values of the dependent functions' Jacobian. It is expected that this matrix be of the correct size (with dimensions n_dependent_variables \(\times\)n_independent_variables).

Definition at line 242 of file ad_drivers.cc.


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