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
Public Types | List of all members
Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType > Class Template Referenceabstract

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

Inheritance diagram for Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >:
[legend]

Public Types

using scalar_type = typename HelperBase< ADNumberTypeCode, ScalarType >::scalar_type
 
using ad_type = typename HelperBase< ADNumberTypeCode, ScalarType >::ad_type
 

Public Member Functions

Constructor / destructor
 CellLevelBase (const unsigned int n_independent_variables, const unsigned int n_dependent_variables)
 
virtual ~CellLevelBase ()=default
 
Independent variables
void register_dof_values (const std::vector< scalar_type > &dof_values)
 
template<typename VectorType >
void register_dof_values (const VectorType &values, const std::vector<::types::global_dof_index > &local_dof_indices)
 
const std::vector< ad_type > & get_sensitive_dof_values () const
 
Operations specific to taped mode: Reusing tapes
void set_dof_values (const std::vector< scalar_type > &dof_values)
 
template<typename VectorType >
void set_dof_values (const VectorType &values, const std::vector<::types::global_dof_index > &local_dof_indices)
 
Dependent variables
virtual void compute_residual (Vector< scalar_type > &residual) const =0
 
virtual void compute_linearization (FullMatrix< scalar_type > &linearization) const =0
 
Interrogation of internal information
std::size_t n_independent_variables () const
 
std::size_t n_dependent_variables () const
 
void print (std::ostream &stream) const
 
void print_values (std::ostream &stream) const
 
void print_tape_stats (const typename Types< ad_type >::tape_index tape_index, std::ostream &stream) const
 
Operations specific to taped mode: Recording tapes
virtual void reset (const unsigned int n_independent_variables=::numbers::invalid_unsigned_int, const unsigned int n_dependent_variables=::numbers::invalid_unsigned_int, const bool clear_registered_tapes=true)
 
bool is_recording () const
 
Types< ad_type >::tape_index active_tape_index () const
 
bool is_registered_tape (const typename Types< ad_type >::tape_index tape_index) const
 
void set_tape_buffer_sizes (const typename Types< ad_type >::tape_buffer_sizes obufsize=64 *1024 *1024, const typename Types< ad_type >::tape_buffer_sizes lbufsize=64 *1024 *1024, const typename Types< ad_type >::tape_buffer_sizes vbufsize=64 *1024 *1024, const typename Types< ad_type >::tape_buffer_sizes tbufsize=64 *1024 *1024)
 
bool start_recording_operations (const typename Types< ad_type >::tape_index tape_index, const bool overwrite_tape=false, const bool keep_independent_values=true)
 
void stop_recording_operations (const bool write_tapes_to_file=false)
 
void activate_recorded_tape (const typename Types< ad_type >::tape_index tape_index)
 
bool recorded_tape_requires_retaping (const typename Types< ad_type >::tape_index tape_index) const
 
bool active_tape_requires_retaping () const
 
void clear_active_tape ()
 

Static Public Member Functions

Operations specific to tapeless mode
static void configure_tapeless_mode (const unsigned int n_independent_variables, const bool ensure_persistent_setting=true)
 

Drivers and taping

TapedDrivers< ad_type, scalar_typetaped_driver
 
TapelessDrivers< ad_type, scalar_typetapeless_driver
 
void activate_tape (const typename Types< ad_type >::tape_index tape_index, const bool read_mode)
 

Independent variables

std::vector< scalar_typeindependent_variable_values
 
std::vector< ad_typeindependent_variables
 
std::vector< boolregistered_independent_variable_values
 
std::vector< boolregistered_marked_independent_variables
 
void reset_registered_independent_variables ()
 
void set_sensitivity_value (const unsigned int index, const scalar_type &value)
 
void mark_independent_variable (const unsigned int index, ad_type &out) const
 
void finalize_sensitive_independent_variables () const
 
void initialize_non_sensitive_independent_variable (const unsigned int index, ad_type &out) const
 
unsigned int n_registered_independent_variables () const
 

Dependent variables

std::vector< ad_typedependent_variables
 
std::vector< boolregistered_marked_dependent_variables
 
void reset_registered_dependent_variables (const bool flag=false)
 
unsigned int n_registered_dependent_variables () const
 
void register_dependent_variable (const unsigned int index, const ad_type &func)
 

Detailed Description

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
class Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >

A general helper class that facilitates the evaluation of a vector of functions, as well as its first derivatives (their Jacobian). This class would typically be used to compute the linearization of a set of local nonlinear equations, but can also be used as the basis of the linearization of the residual vector defined on the level of a finite element (for example, in order to compute the Jacobian matrix necessary in Newton-type solvers for nonlinear problems).

Note
When using the cell-level taped AD methods in 3d and/or with higher order elements, it is incredibly easy to exceed the tape buffer size. The reason for this is two-fold:
  1. there are many independent variables (the local degrees-of-freedom) to take the derivatives with respect to, and
  2. the expressions for the dependent variables (each being a component of the residual vector) in terms of all of the independent variables are lengthy, especially when non-trivial constitutive laws are considered. These buffer variables dictate the amount of memory allocated to a tape before it is written to file (at a significant performance loss). Therefore for ADOL-C taped AD numbers, it may be desirable to create a file ".adolcrc" in the program run directory and set the buffer size therein (as is suggested by the ADOL-C manual). For example, the following settings increase the default buffer size by 128 times:
    "OBUFSIZE" "67108864"
    "LBUFSIZE" "67108864"
    "VBUFSIZE" "67108864"
    "TBUFSIZE" "67108864"
    Note that the quotation marks are mandatory. An alternative approach that allows for run-time decision making is to use the HelperBase::set_tape_buffer_sizes() function before starting taping (as done via the HelperBase::start_recording_operations() function).
Warning
ADOL-C does not support the standard threading models used by deal.II, so this class should not be embedded within a multithreaded function when using ADOL-C number types. It is, however, suitable for use in both serial and MPI routines.

Definition at line 835 of file ad_helpers.h.

Member Typedef Documentation

◆ scalar_type

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
using Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >::scalar_type = typename HelperBase<ADNumberTypeCode, ScalarType>::scalar_type

Type definition for the floating point number type that is used in, and results from, all computations.

Definition at line 842 of file ad_helpers.h.

◆ ad_type

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
using Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >::ad_type = typename HelperBase<ADNumberTypeCode, ScalarType>::ad_type

Type definition for the auto-differentiation number type that is used in all computations.

Definition at line 849 of file ad_helpers.h.

Constructor & Destructor Documentation

◆ CellLevelBase()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >::CellLevelBase ( const unsigned int  n_independent_variables,
const unsigned int  n_dependent_variables 
)

The constructor for the class.

Parameters
[in]n_independent_variablesThe number of independent variables that will be used in the definition of the functions that it is desired to compute the sensitivities of. In the computation of \(\mathbf{f}(\mathbf{X})\), this will be the number of inputs \(\mathbf{X}\), i.e., the dimension of the domain space.
[in]n_dependent_variablesThe number of scalar functions to be defined that will have a sensitivity to the given independent variables. In the computation of \(\mathbf{f}(\mathbf{X})\), this will be the number of outputs \(\mathbf{f}\), i.e., the dimension of the image space.

Definition at line 710 of file ad_helpers.cc.

◆ ~CellLevelBase()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
virtual Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >::~CellLevelBase ( )
virtualdefault

Destructor

Member Function Documentation

◆ register_dof_values() [1/2]

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >::register_dof_values ( const std::vector< scalar_type > &  dof_values)

Register the complete set of independent variables \(\mathbf{X}\) that represent the local degree of freedom values.

Parameters
[in]dof_valuesA vector field associated with local degree of freedom values on the current finite element. These define the values of all independent variables. When considering taped AD numbers with branching functions, to avoid potential issues with branch switching it may be a good idea to choose these values close or equal to those that will be later evaluated and linearized around.
Note
The input value type must correspond to this class's scalar_type. Depending on the selected ADNumberTypeCode, this may or may not correspond with the ScalarType prescribed as a template argument.
For taped AD numbers, this operation is only valid in recording mode.

Definition at line 721 of file ad_helpers.cc.

◆ register_dof_values() [2/2]

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
template<typename VectorType >
void Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >::register_dof_values ( const VectorType &  values,
const std::vector<::types::global_dof_index > &  local_dof_indices 
)

Register the complete set of independent variables \(\mathbf{X}\) that represent the local degree of freedom values.

Parameters
[in]valuesA global field from which the values of all independent variables will be extracted. This typically will be the solution vector around which point a residual vector is to be computed and around which linearization is to occur. When considering taped AD numbers with branching functions, to avoid potential issues with branch switching it may be a good idea to choose these values close or equal to those that will be later evaluated and linearized around.
[in]local_dof_indicesA vector of degree of freedom indices from which to extract the local degree of freedom values. This would typically obtained by calling cell->get_dof_indices().
Note
For taped AD numbers, this operation is only valid in recording mode.

◆ get_sensitive_dof_values()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
const std::vector< typename CellLevelBase< ADNumberTypeCode, ScalarType >::ad_type > & Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >::get_sensitive_dof_values

Return the complete set of degree of freedom values as represented by auto-differentiable numbers. These are the independent variables \(\mathbf{X}\) about which the solution is linearized.

This function indicates to the AD library that implements the auto-differentiable number type that operations performed on these numbers are to be tracked so they are considered "sensitive" variables. This is, therefore, the set of variables with which one would then perform computations, and based on which one can then extract both the value of the function and its derivatives with the member functions below. The values of the components of the returned object are initialized to the values set with register_independent_variable().

Returns
An array of auto-differentiable type numbers representing the local degree of freedom values.
Note
For taped AD numbers, this operation is only valid in recording mode.

Definition at line 744 of file ad_helpers.cc.

◆ set_dof_values() [1/2]

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >::set_dof_values ( const std::vector< scalar_type > &  dof_values)

Set the values for the independent variables \(\mathbf{X}\), i.e., the linearization point.

Parameters
[in]dof_valuesA vector field associated with local degree of freedom values on the current finite element. These define the values of all independent variables.
Note
The input value type must correspond to this class's scalar_type. Depending on the selected ADNumberTypeCode, this may or may not correspond with the ScalarType prescribed as a template argument.
If the keep_independent_values flag has been set when HelperBase::start_recording_operations() is called then the tape is immediately usable after creation, and the values of the independent variables set by register_dof_values() are those at which the function is to be evaluated. In this case, a separate call to this function is not strictly necessary.

Definition at line 769 of file ad_helpers.cc.

◆ set_dof_values() [2/2]

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
template<typename VectorType >
void Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >::set_dof_values ( const VectorType &  values,
const std::vector<::types::global_dof_index > &  local_dof_indices 
)

Set the values for the independent variables \(\mathbf{X}\), i.e., the linearization point.

Parameters
[in]valuesA vector field from which the values of all independent variables is to be extracted.
[in]local_dof_indicesA vector of degree of freedom indices from which to extract the local degree of freedom values. This would typically obtained by calling cell->get_dof_indices().
Note
If the keep_independent_values flag has been set when HelperBase::start_recording_operations() is called then the tape is immediately usable after creation, and the values of the independent variables set by register_dof_values() are those at which the function is to be evaluated. In this case, a separate call to this function is not strictly necessary.

◆ compute_residual()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
virtual void Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >::compute_residual ( Vector< scalar_type > &  residual) const
pure virtual

Compute the value of the residual vector field \(\mathbf{r}(\mathbf{X})\).

Parameters
[out]residualA Vector object with the value for each component of the vector field evaluated at the point defined by the independent variable values.
Note
The size of the residual vector is determined by the derived classes, as it depends on the order of the dependent variable(s) derivative(s) that it represents. Code examples that show how to use this interface will be provided in the documentation of the derived classes.

Implemented in Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >, and Differentiation::AD::ResidualLinearization< ADNumberTypeCode, ScalarType >.

◆ compute_linearization()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
virtual void Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >::compute_linearization ( FullMatrix< scalar_type > &  linearization) const
pure virtual

Compute the gradient (first derivative) of the residual vector field with respect to all independent variables, i.e.

\[ \frac{\partial\mathbf{r}(\mathbf{X})}{\partial\mathbf{X}} \]

Parameters
[out]linearizationA FullMatrix with the gradient of each component of the vector field evaluated at the point defined by the independent variable values.
Note
The dimensions of the linearization matrix is determined by the derived classes, as it depends on the order of the dependent variable(s) derivative(s) that it represents. Code examples that show how to use this interface will be provided in the documentation of the derived classes.

Implemented in Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >, and Differentiation::AD::ResidualLinearization< ADNumberTypeCode, ScalarType >.

◆ n_independent_variables()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
std::size_t Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::n_independent_variables
inherited

Return the number of independent variables that this object expects to work with. This is the dimension of the domain space.

Definition at line 241 of file ad_helpers.cc.

◆ n_dependent_variables()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
std::size_t Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::n_dependent_variables
inherited

Return the number of dependent variables that this object expects to operate on. This is the dimension of the image space.

Definition at line 262 of file ad_helpers.cc.

◆ print()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::print ( std::ostream &  stream) const
inherited

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 309 of file ad_helpers.cc.

◆ print_values()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::print_values ( std::ostream &  stream) const
inherited

Print the values currently assigned to the independent variables.

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

Definition at line 359 of file ad_helpers.cc.

◆ print_tape_stats()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::print_tape_stats ( const typename Types< ad_type >::tape_index  tape_index,
std::ostream &  stream 
) const
inherited

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.
Note
This function only produces meaningful output when ad_type is a taped auto-differentiable number.

Definition at line 373 of file ad_helpers.cc.

◆ configure_tapeless_mode()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::configure_tapeless_mode ( const unsigned int  n_independent_variables,
const bool  ensure_persistent_setting = true 
)
staticinherited

Pre-specify the number of independent_variables to be used in tapeless mode.

Although this function is called internally in the HelperBase constructor, there may be occasions when ADOL-C tapeless numbers (adtl::adoubles) are created before an instance of this class is created. This function therefore allows one to declare at the earliest possible instance how many directional derivatives will be considered in tapeless mode.

Warning
With ensure_persistent_setting set to true when the ad_type is an ADOL-C tapeless number, calling this function leaves the set number of directional derivatives in a persistent state. It will therefore not be possible to further modify the number of directional derivatives to be tracked by adtl::adoubles's during course of the program's execution.

Definition at line 448 of file ad_helpers.cc.

◆ reset()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::reset ( const unsigned int  n_independent_variables = ::numbers::invalid_unsigned_int,
const unsigned int  n_dependent_variables = ::numbers::invalid_unsigned_int,
const bool  clear_registered_tapes = true 
)
virtualinherited

Reset the state of the helper class.

When an instance of an HelperBase is stored as a class member object with the intention to reuse its instance, it may be necessary to reset the state of the object before use. This is because, internally, there is error checking performed to ensure that the correct auto-differentiable data is being tracked and used only when appropriate. This function clears all member data and, therefore, allows the state of all internal flags to be safely reset to their initial state.

In the rare case that the number of independent or dependent variables has changed, this can also reconfigured by passing in the appropriate arguments to the function.

Parameters
[in]n_independent_variablesThe number of independent variables that will be used in the definition of the functions that it is desired to compute the sensitivities of. In the computation of \(\mathbf{f}(\mathbf{X})\), this will be the number of inputs \(\mathbf{X}\), i.e., the dimension of the domain space.
[in]n_dependent_variablesThe number of scalar functions to be defined that will have a sensitivity to the given independent variables. In the computation of \(\mathbf{f}(\mathbf{X})\), this will be the number of outputs \(\mathbf{f}\), i.e., the dimension of the image space.
[in]clear_registered_tapesA flag that indicates the that list of registered_tapes must be cleared. If set to true then the data structure that tracks which tapes have been recorded is cleared as well. It is then expected that any preexisting tapes be re-recorded.
Note
This also resets the active tape number to an invalid number, and deactivates the recording mode for taped variables.

Reimplemented in Differentiation::AD::PointLevelFunctionsBase< dim, ADNumberTypeCode, ScalarType >, and Differentiation::AD::PointLevelFunctionsBase< dim, ADNumberTypeCode, double >.

Definition at line 390 of file ad_helpers.cc.

◆ is_recording()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
bool Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::is_recording
inherited

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

Definition at line 271 of file ad_helpers.cc.

◆ active_tape_index()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
Types< typenameHelperBase< ADNumberTypeCode, ScalarType >::ad_type >::tape_index Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::active_tape_index
inherited

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

Definition at line 284 of file ad_helpers.cc.

◆ is_registered_tape()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
bool Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::is_registered_tape ( const typename Types< ad_type >::tape_index  tape_index) const
inherited

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

Definition at line 296 of file ad_helpers.cc.

◆ set_tape_buffer_sizes()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::set_tape_buffer_sizes ( const typename Types< ad_type >::tape_buffer_sizes  obufsize = 64 * 1024 * 1024,
const typename Types< ad_type >::tape_buffer_sizes  lbufsize = 64 * 1024 * 1024,
const typename Types< ad_type >::tape_buffer_sizes  vbufsize = 64 * 1024 * 1024,
const typename Types< ad_type >::tape_buffer_sizes  tbufsize = 64 * 1024 * 1024 
)
inherited

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 559 of file ad_helpers.cc.

◆ start_recording_operations()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
bool Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::start_recording_operations ( const typename Types< ad_type >::tape_index  tape_index,
const bool  overwrite_tape = false,
const bool  keep_independent_values = true 
)
inherited

Enable recording mode for a given tape. The use of this function is mandatory if the auto-differentiable number is a taped type. However, for the purpose of developing generic code, it can also be safely called for tapeless auto-differentiable numbers.

The operations that take place between this function call and that of stop_recording_operations() are recorded to the tape and can be replayed and reevaluated as necessary.

The typical set of operations to be performed during this "recording" phase (between the calls to start_recording_operations() and stop_recording_operations() ) are:

  • Definition of some independent variables via register_independent_variable() or register_independent_variables(). These define the branch of operations tracked by the tape. If the keep flag is set to true then these represent precisely the point about which the function derivatives are to be computed. If the keep flag is set to false then these only represent dummy values, and the point at which the function derivatives are to be computed must be set by calling set_independent_variables() again.
  • Extraction of a set of independent variables of auto-differentiable type using get_sensitive_variables(). These are then tracked during later computations.
  • Defining the dependent variables via register_dependent_variable() or register_dependent_variables(). These are the functions that will be differentiated with respect to the independent variables.
Parameters
[in]tape_indexThe index of the tape to be written
[in]overwrite_tapeExpress whether tapes are allowed to be overwritten. If true then any existing tape with a given tape_index will be destroyed and a new tape traced over it.
[in]keep_independent_valuesDetermines whether the numerical values of all independent variables are recorded in the tape buffer. If true, then the tape can be immediately used to perform computations after recording is complete.
Note
During the recording phase, no value(), gradient(), hessian(), or jacobian() operations can be performed.
The chosen tape index must be greater than Numbers<ad_type>::invalid_tape_index and less than Numbers<ad_type>::max_tape_index.

Definition at line 578 of file ad_helpers.cc.

◆ stop_recording_operations()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::stop_recording_operations ( const bool  write_tapes_to_file = false)
inherited

Disable recording mode for a given tape. The use of this function is mandatory if the auto-differentiable number is a taped type. However, for the purpose of developing generic code, it can also be safely called for tapeless auto-differentiable numbers.

Note
After this function call, the tape is considered ready for use and operations such as value(), gradient() or hessian() can be executed.
For taped AD numbers, this operation is only valid in recording mode.

Definition at line 643 of file ad_helpers.cc.

◆ activate_recorded_tape()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::activate_recorded_tape ( const typename Types< ad_type >::tape_index  tape_index)
inherited

Select a pre-recorded tape to read from.

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

Definition at line 474 of file ad_helpers.cc.

◆ recorded_tape_requires_retaping()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
bool Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::recorded_tape_requires_retaping ( const typename Types< ad_type >::tape_index  tape_index) const
inherited

Return a flag that, when true, indicates that the retaping of the dependent function is necessary for a reliable computation to be performed on a tape with the given tape_index. 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.

For the output of this function to be meaningful, it must be called after activate_recorded_tape() is called and the new evaluation point for the tape (i.e. values of the independent variables) have been set and subsequently used (i.e. in the determination of the values or derivatives of the dependent variables).

Definition at line 484 of file ad_helpers.cc.

◆ active_tape_requires_retaping()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
bool Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::active_tape_requires_retaping
inherited

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.

For the output of this function to be meaningful, it must be called after activate_recorded_tape() is called and the new evaluation point for the tape (i.e. values of the independent variables) have been set and subsequently used (i.e. in the determination of the values or derivatives of the dependent variables).

Definition at line 497 of file ad_helpers.cc.

◆ clear_active_tape()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::clear_active_tape
inherited

Clears and removes the currently active tape.

This is typically only necessary when branch switching is detected on the original tape at evaluation point. This state can be checked using the active_tape_requires_retaping() function.

Definition at line 510 of file ad_helpers.cc.

◆ activate_tape()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::activate_tape ( const typename Types< ad_type >::tape_index  tape_index,
const bool  read_mode 
)
protectedinherited

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.
[in]read_modeA flag that marks whether or not we expect to read data from a preexisting tape.
Note
The chosen tape index must be greater than Numbers<ad_type>::invalid_tape_index and less than Numbers<ad_type>::max_tape_index.

Definition at line 522 of file ad_helpers.cc.

◆ reset_registered_independent_variables()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::reset_registered_independent_variables
protectedinherited

Reset the boolean vector registered_independent_variable_values that indicates which independent variables we've been manipulating for the current set of operations.

Definition at line 82 of file ad_helpers.cc.

◆ set_sensitivity_value()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::set_sensitivity_value ( const unsigned int  index,
const scalar_type value 
)
protectedinherited

Set the actual value of the independent variable \(X_{i}\).

Parameters
[in]indexThe index in the vector of independent variables.
[in]valueThe value to set the index'd independent variable to.

Definition at line 105 of file ad_helpers.cc.

◆ mark_independent_variable()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::mark_independent_variable ( const unsigned int  index,
ad_type out 
) const
protectedinherited

Initialize an independent variable \(X_{i}\) such that subsequent operations performed with it are tracked.

Note
Care must be taken to mark each independent variable only once.
The order in which the independent variables are marked defines the order of all future internal operations. They must be manipulated in the same order as that in which they are first marked. If not then, for example, ADOL-C won't throw an error, but rather it might complain nonsensically during later computations or produce garbage results.
Parameters
[in]indexThe index in the vector of independent variables.
[out]outAn auto-differentiable number that is ready for use in computations. The operations that are performed with it are recorded on the tape and will be considered in the computation of dependent variable values.

Definition at line 142 of file ad_helpers.cc.

◆ finalize_sensitive_independent_variables()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::finalize_sensitive_independent_variables
protectedinherited

Finalize the state of the independent variables before use.

This step and the storage of the independent variables is done separately because some derived classes may offer the capability to add independent variables in a staggered manner. This function is to be triggered when these values are considered finalized and we can safely initialize the sensitive equivalents of those values.

Definition at line 181 of file ad_helpers.cc.

◆ initialize_non_sensitive_independent_variable()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::initialize_non_sensitive_independent_variable ( const unsigned int  index,
ad_type out 
) const
protectedinherited

Initialize an independent variable \(X_{i}\).

Parameters
[out]outAn auto-differentiable number that is ready for use in standard computations. The operations that are performed with it are not recorded on the tape, and so should only be used when not in recording mode.
[in]indexThe index in the vector of independent variables.

Definition at line 204 of file ad_helpers.cc.

◆ n_registered_independent_variables()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
unsigned int Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::n_registered_independent_variables
protectedinherited

The number of independent variables that have been manipulated within a set of operations.

Definition at line 230 of file ad_helpers.cc.

◆ reset_registered_dependent_variables()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::reset_registered_dependent_variables ( const bool  flag = false)
protectedinherited

Reset the boolean vector registered_marked_dependent_variables that indicates which independent variables have been manipulated by the current set of operations. All entries in the vector are set to the value of the flag.

Definition at line 93 of file ad_helpers.cc.

◆ n_registered_dependent_variables()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
unsigned int Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::n_registered_dependent_variables
protectedinherited

The number of dependent variables that have been registered.

Definition at line 250 of file ad_helpers.cc.

◆ register_dependent_variable()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::register_dependent_variable ( const unsigned int  index,
const ad_type func 
)
protectedinherited

Register the definition of the index'th dependent variable \(f(\mathbf{X})\).

Parameters
[in]indexThe index of the entry in the global list of dependent variables that this function belongs to.
[in]funcThe recorded function that defines a dependent variable.
Note
Each dependent variable must only be registered once.

Definition at line 679 of file ad_helpers.cc.

Member Data Documentation

◆ taped_driver

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
TapedDrivers<ad_type, scalar_type> Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::taped_driver
protectedinherited

An object used to help manage stored tapes.

In the event that the ad_type is a tapeless AD type, then the object constructed here is, effectively, a dummy one.

Definition at line 591 of file ad_helpers.h.

◆ tapeless_driver

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
TapelessDrivers<ad_type, scalar_type> Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::tapeless_driver
protectedinherited

An object used to help manage tapeless data structures.

In the event that the ad_type is a taped AD type, then the object constructed here is, effectively, a dummy one.

Definition at line 599 of file ad_helpers.h.

◆ independent_variable_values

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
std::vector<scalar_type> Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::independent_variable_values
mutableprotectedinherited

A set of independent variables \(\mathbf{X}\) that differentiation will be performed with respect to.

The gradients and Hessians of dependent variables will be computed at these finite values.

Definition at line 635 of file ad_helpers.h.

◆ independent_variables

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
std::vector<ad_type> Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::independent_variables
mutableprotectedinherited

A set of sensitive independent variables \(\mathbf{X}\) that differentiation will be performed with respect to.

The gradients and Hessians of dependent variables will be computed using these configured AD numbers. Note that only reverse-mode AD requires that the sensitive independent variables be stored.

Definition at line 645 of file ad_helpers.h.

◆ registered_independent_variable_values

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
std::vector<bool> Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::registered_independent_variable_values
protectedinherited

A list of registered independent variables that have been manipulated for a given set of operations.

Definition at line 651 of file ad_helpers.h.

◆ registered_marked_independent_variables

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
std::vector<bool> Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::registered_marked_independent_variables
mutableprotectedinherited

A list of registered independent variables that have been extracted and their sensitivities marked.

Definition at line 657 of file ad_helpers.h.

◆ dependent_variables

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
std::vector<ad_type> Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::dependent_variables
protectedinherited

The set of dependent variables \(\mathbf{f}(\mathbf{X})\) of which the derivatives with respect to \(\mathbf{X}\) will be computed.

The gradients and Hessians of these dependent variables will be computed at the values \(\mathbf{X}\) set with the set_sensitivity_value() function.

Note
These are stored as an ad_type so that we can use them to compute function values and directional derivatives in the case that tapeless numbers are used

Definition at line 749 of file ad_helpers.h.

◆ registered_marked_dependent_variables

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
std::vector<bool> Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >::registered_marked_dependent_variables
protectedinherited

A list of registered dependent variables.

Definition at line 754 of file ad_helpers.h.


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