deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
|
#include <deal.II/differentiation/ad/ad_helpers.h>
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 | |
ScalarFunction (const unsigned int n_independent_variables) | |
virtual | ~ScalarFunction ()=default |
Operations specific to taped mode: Reusing tapes | |
void | set_independent_variables (const std::vector< scalar_type > &values) |
template<typename ValueType , typename ExtractorType > | |
void | set_independent_variable (const ValueType &value, const ExtractorType &extractor) |
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 | |
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) |
Static Public Attributes | |
static constexpr unsigned int | dimension = dim |
Dependent variables | |
void | register_dependent_variable (const ad_type &func) |
scalar_type | compute_value () const |
void | compute_gradient (Vector< scalar_type > &gradient) const |
void | compute_hessian (FullMatrix< scalar_type > &hessian) const |
template<typename ExtractorType_Row > | |
static internal::ScalarFieldGradient< dim, scalar_type, ExtractorType_Row >::type | extract_gradient_component (const Vector< scalar_type > &gradient, const ExtractorType_Row &extractor_row) |
template<typename ExtractorType_Row , typename ExtractorType_Col > | |
static internal::ScalarFieldHessian< dim, scalar_type, ExtractorType_Row, ExtractorType_Col >::type | extract_hessian_component (const FullMatrix< scalar_type > &hessian, const ExtractorType_Row &extractor_row, const ExtractorType_Col &extractor_col) |
static Tensor< 0, dim, scalar_type > | extract_hessian_component (const FullMatrix< scalar_type > &hessian, const FEValuesExtractors::Scalar &extractor_row, const FEValuesExtractors::Scalar &extractor_col) |
static SymmetricTensor< 4, dim, scalar_type > | extract_hessian_component (const FullMatrix< scalar_type > &hessian, const FEValuesExtractors::SymmetricTensor< 2 > &extractor_row, const FEValuesExtractors::SymmetricTensor< 2 > &extractor_col) |
Independent variables | |
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) override |
void | register_independent_variables (const std::vector< scalar_type > &values) |
template<typename ValueType , typename ExtractorType > | |
void | register_independent_variable (const ValueType &value, const ExtractorType &extractor) |
const std::vector< ad_type > & | get_sensitive_variables () const |
template<typename ExtractorType > | |
internal::Extractor< dim, ExtractorType >::template tensor_type< typename HelperBase< ADNumberTypeCode, ScalarType >::ad_type > | get_sensitive_variables (const ExtractorType &extractor) const |
std::vector< bool > | symmetric_independent_variables |
void | set_sensitivity_value (const unsigned int index, const bool symmetric_component, const scalar_type &value) |
bool | is_symmetric_independent_variable (const unsigned int index) const |
unsigned int | n_symmetric_independent_variables () const |
Independent variables | |
void | set_sensitivity_value (const unsigned int index, const scalar_type &value) |
void | reset_registered_independent_variables () |
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 |
std::vector< scalar_type > | independent_variable_values |
std::vector< ad_type > | independent_variables |
std::vector< bool > | registered_independent_variable_values |
std::vector< bool > | registered_marked_independent_variables |
Drivers and taping | |
void | activate_tape (const typename Types< ad_type >::tape_index tape_index, const bool read_mode) |
TapedDrivers< ad_type, scalar_type > | taped_driver |
TapelessDrivers< ad_type, scalar_type > | tapeless_driver |
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) |
std::vector< ad_type > | dependent_variables |
std::vector< bool > | registered_marked_dependent_variables |
A helper class that facilitates the evaluation of a scalar function, its first derivatives (gradient), and its second derivatives (Hessian). This class would typically be used to compute the first and second derivatives of a stored energy function defined at a quadrature point. It can also be used to compute derivatives of any other scalar field so long as all its dependencies on the independent variables are explicit (that is to say, no independent variables may have some implicit dependence on one another).
An example of its usage in the case of a multi-field constitutive law might be as follows:
Definition at line 3116 of file ad_helpers.h.
using Differentiation::AD::ScalarFunction< dim, 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 3124 of file ad_helpers.h.
using Differentiation::AD::ScalarFunction< dim, 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 3131 of file ad_helpers.h.
Differentiation::AD::ScalarFunction< dim, ADNumberTypeCode, ScalarType >::ScalarFunction | ( | const unsigned int | n_independent_variables | ) |
The constructor for the class.
[in] | n_independent_variables | The 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. |
Definition at line 1320 of file ad_helpers.cc.
|
virtualdefault |
Destructor.
void Differentiation::AD::ScalarFunction< dim, ADNumberTypeCode, ScalarType >::register_dependent_variable | ( | const ad_type & | func | ) |
Register the definition of the scalar field \(\Psi(\mathbf{X})\).
[in] | func | The recorded function that defines a dependent variable. |
Definition at line 1333 of file ad_helpers.cc.
ScalarFunction< dim, ADNumberTypeCode, ScalarType >::scalar_type Differentiation::AD::ScalarFunction< dim, ADNumberTypeCode, ScalarType >::compute_value | ( | ) | const |
Compute the value of the scalar field \(\Psi(\mathbf{X})\) using the tape as opposed to executing the source code.
Definition at line 1347 of file ad_helpers.cc.
void Differentiation::AD::ScalarFunction< dim, ADNumberTypeCode, ScalarType >::compute_gradient | ( | Vector< scalar_type > & | gradient | ) | const |
Compute the gradient (first derivative) of the scalar field with respect to all independent variables, i.e.
\[ \frac{\partial\Psi(\mathbf{X})}{\partial\mathbf{X}} \]
[out] | gradient | A Vector with the values for the scalar field gradient (first derivatives) evaluated at the point defined by the independent variable values. The output gradient vector has a length corresponding to n_independent_variables . |
Definition at line 1397 of file ad_helpers.cc.
void Differentiation::AD::ScalarFunction< dim, ADNumberTypeCode, ScalarType >::compute_hessian | ( | FullMatrix< scalar_type > & | hessian | ) | const |
Compute the Hessian (second derivative) of the scalar field with respect to all independent variables, i.e.
\[ \frac{\partial^{2}\Psi(\mathbf{X})}{\partial\mathbf{X} \otimes \partial\mathbf{X}} \]
[out] | hessian | A FullMatrix with the values for the scalar field Hessian (second derivatives) evaluated at the point defined by the independent variable values. The output hessian matrix has dimensions corresponding to n_independent_variables \(\times\)n_independent_variables . |
Definition at line 1471 of file ad_helpers.cc.
|
static |
Extract the function gradient for a subset of independent variables \(\mathbf{A} \subset \mathbf{X}\), i.e.
\[ \frac{\partial\Psi(\mathbf{X})}{\partial\mathbf{A}} \]
[in] | gradient | The gradient of the scalar function with respect to all independent variables, i.e., that returned by compute_gradient(). |
[in] | extractor_row | An extractor associated with the input field variables. This effectively defines which components of the global set of independent variables this field is associated with. |
extractor_row
. This corresponds to subsetting a whole set of rows of the gradient vector, scaling those entries to take account of tensor symmetries, and then reshaping the (sub-)vector so obtained into a tensor, the final result. For example, if extractor_row
is a FEValuesExtractors::Vector and extractor_col
is a FEValuesExtractors::Tensor, then the returned object is a Tensor of rank 3, with its first index associated with the field corresponding to the row extractor and the second and third indices associated with the field corresponding to the column extractor. Similarly, if extractor_row
is a FEValuesExtractors::SymmetricTensor and extractor_col
is a FEValuesExtractors::SymmetricTensor, then the returned object is a SymmetricTensor of rank 4, with its first two indices associated with the field corresponding to the row extractor and the last two indices associated with the field corresponding to the column extractor.
|
static |
Extract the function Hessian for a subset of independent variables \(\mathbf{A},\mathbf{B} \subset \mathbf{X}\), i.e.
\[ \frac{}{\partial\mathbf{B}} \left[ \frac{\partial\Psi(\mathbf{X})}{\partial\mathbf{A}} \right] = \frac{\partial^{2}\Psi(\mathbf{X})}{\partial\mathbf{B} \otimes \partial\mathbf{A}} \]
[in] | hessian | The Hessian of the scalar function with respect to all independent variables, i.e., that returned by compute_hessian(). |
[in] | extractor_row | An extractor associated with the input field variables for which the first index of the Hessian is extracted. |
[in] | extractor_col | An extractor associated with the input field variables for which the second index of the Hessian is extracted. |
extractor_row
and extractor_col
. This corresponds to subsetting a whole set of rows and columns of the Hessian matrix, scaling those entries to take account of tensor symmetries, and then reshaping the (sub-)matrix so obtained into a tensor, the final result. For example, if extractor_row
is a FEValuesExtractors::Vector and extractor_col
is a FEValuesExtractors::Tensor, then the returned object is a Tensor of rank 3, with its first index associated with the field corresponding to the row extractor and the second and third indices associated with the field corresponding to the column extractor. Similarly, if extractor_row
is a FEValuesExtractors::SymmetricTensor and extractor_col
is a FEValuesExtractors::SymmetricTensor, then the returned object is a SymmetricTensor of rank 4, with its first two indices associated with the field corresponding to the row extractor and the last two indices associated with the field corresponding to the column extractor.
|
static |
Extract the function Hessian for a subset of independent variables \(\mathbf{A},\mathbf{B} \subset \mathbf{X}\), i.e.
\[ \frac{}{\partial\mathbf{B}} \left[ \frac{\partial\Psi(\mathbf{X})}{\partial\mathbf{A}} \right] \]
This function is a specialization of the above for rank-0 tensors (scalars). This corresponds to extracting a single entry of the Hessian matrix because both extractors imply selection of just a single row or column of the matrix.
Definition at line 1569 of file ad_helpers.cc.
|
static |
Extract the function Hessian for a subset of independent variables \(\mathbf{A},\mathbf{B} \subset \mathbf{X}\), i.e.
\[ \frac{}{\partial\mathbf{B}} \left[ \frac{\partial\Psi(\mathbf{X})}{\partial\mathbf{A}} \right] \]
This function is a specialization of the above for rank-4 symmetric tensors.
Definition at line 1606 of file ad_helpers.cc.
|
overridevirtualinherited |
Reset the state of the class.
When an instance of the class 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 be reconfigured by passing in the appropriate arguments to the function.
[in] | n_independent_variables | The 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_variables | The 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_tapes | A 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. |
Reimplemented from Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >.
Definition at line 1167 of file ad_helpers.cc.
|
inherited |
Register the complete set of independent variables \(\mathbf{X}\).
[in] | values | A field that defines 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 differentiated around. |
scalar_type
. Depending on the selected ADNumberTypeCode
, this may or may not correspond with the ScalarType
prescribed as a template argument.Definition at line 1218 of file ad_helpers.cc.
|
inherited |
Register the subset of independent variables \(\mathbf{A} \subset \mathbf{X}\).
[in] | value | A field that defines a number of 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 differentiated around. |
[in] | extractor | An extractor associated with the input field variables. This effectively defines which components of the global set of independent variables this field is associated with. |
scalar_type
. Depending on the selected ADNumberTypeCode
, this may or may not correspond with the ScalarType
prescribed as a template argument.ValueType
. So, for example, if a value is a rank-1 tensor (i.e. of type Tensor<1,dim,scalar_type>), then the extractor must be an FEValuesExtractors::Vector or FEValuesExtractors::Tensor<1>.
|
inherited |
Return the complete set of independent variables as represented by auto-differentiable numbers. These are the independent variables \(\mathbf{X}\) at which the dependent values are evaluated and differentiated.
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().
Definition at line 1244 of file ad_helpers.cc.
|
inherited |
|
inherited |
Set the values for the independent variables \(\mathbf{X}\).
[in] | values | A vector that defines the values of all independent variables. |
scalar_type
. Depending on the selected ADNumberTypeCode
, this may or may not correspond with the ScalarType
prescribed as a template argument.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_independent_variables() 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 1294 of file ad_helpers.cc.
|
inherited |
Set the values for a subset of independent variables \(\mathbf{A} \subset \mathbf{X}\).
[in] | value | A field that defines the values of a number of independent variables. |
[in] | extractor | An extractor associated with the input field variables. This effectively defines which components of the global set of independent variables this field is associated with. |
scalar_type
. Depending on the selected ADNumberTypeCode
, this may or may not correspond with the ScalarType
prescribed as a template argument.ValueType
. So, for example, if a value is a rank-1 tensor (i.e. of type Tensor<1,dim,scalar_type>), then the extractor must be an FEValuesExtractors::Vector or FEValuesExtractors::Tensor<1>.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_independent_variable() are those at which the function is to be evaluated. In this case, a separate call to this function is not strictly necessary.
|
protectedinherited |
Set the actual value of the independent variable \(X_{i}\).
[in] | index | The index in the vector of independent variables. |
[in] | symmetric_component | Mark whether this index relates to a component of a field that has a symmetric counterpart (e.g. if index represents an off-diagonal entry in a symmetric tensor). |
[in] | value | The value to set the index'd independent variable to. |
Definition at line 1272 of file ad_helpers.cc.
|
protectedinherited |
Set the actual value of the independent variable \(X_{i}\).
[in] | index | The index in the vector of independent variables. |
[in] | value | The value to set the index'd independent variable to. |
Definition at line 104 of file ad_helpers.cc.
|
protectedinherited |
Return whether the index'th
independent variables is one for which we must take into account symmetry when extracting their gradient or Hessian values.
Definition at line 1190 of file ad_helpers.cc.
|
protectedinherited |
Return the number of independent variables that have been marked as being components of a symmetric field.
Definition at line 1204 of file ad_helpers.cc.
|
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 240 of file ad_helpers.cc.
|
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 261 of file ad_helpers.cc.
|
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.
[in] | stream | The output stream to which the values are to be written. |
Definition at line 308 of file ad_helpers.cc.
|
inherited |
Print the values currently assigned to the independent variables.
[in] | stream | The output stream to which the values are to be written. |
Definition at line 358 of file ad_helpers.cc.
|
inherited |
Print the statistics regarding the usage of the tapes.
[in] | tape_index | The index of the tape to get the statistics of. |
[out] | stream | The output stream to which the values are to be written. |
ad_type
is a taped auto-differentiable number. Definition at line 372 of file ad_helpers.cc.
|
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.
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 447 of file ad_helpers.cc.
|
inherited |
Return whether or not this class is tracking calculations performed with its marked independent variables.
Definition at line 270 of file ad_helpers.cc.
|
inherited |
Return the tape index which is currently activated for recording or reading.
Definition at line 283 of file ad_helpers.cc.
|
inherited |
Return whether or not a tape number has already been used or registered.
Definition at line 295 of file ad_helpers.cc.
|
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.
[in] | obufsize | ADOL-C operations buffer size |
[in] | lbufsize | ADOL-C locations buffer size |
[in] | vbufsize | ADOL-C value buffer size |
[in] | tbufsize | ADOL-C Taylor buffer size |
Definition at line 558 of file ad_helpers.cc.
|
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:
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.[in] | tape_index | The index of the tape to be written |
[in] | overwrite_tape | Express 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_values | Determines 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. |
Definition at line 577 of file ad_helpers.cc.
|
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.
Definition at line 642 of file ad_helpers.cc.
|
inherited |
Select a pre-recorded tape to read from.
[in] | tape_index | The index of the tape to be read from. |
Definition at line 473 of file ad_helpers.cc.
|
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:
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 483 of file ad_helpers.cc.
|
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:
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 496 of file ad_helpers.cc.
|
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 509 of file ad_helpers.cc.
|
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.
[in] | tape_index | The index of the tape to be written to/read from. |
[in] | read_mode | A flag that marks whether or not we expect to read data from a preexisting tape. |
Definition at line 521 of file ad_helpers.cc.
|
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 81 of file ad_helpers.cc.
|
protectedinherited |
Initialize an independent variable \(X_{i}\) such that subsequent operations performed with it are tracked.
[in] | index | The index in the vector of independent variables. |
[out] | out | An 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 141 of file ad_helpers.cc.
|
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 180 of file ad_helpers.cc.
|
protectedinherited |
Initialize an independent variable \(X_{i}\).
[out] | out | An 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] | index | The index in the vector of independent variables. |
Definition at line 203 of file ad_helpers.cc.
|
protectedinherited |
The number of independent variables that have been manipulated within a set of operations.
Definition at line 229 of file ad_helpers.cc.
|
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 92 of file ad_helpers.cc.
|
protectedinherited |
The number of dependent variables that have been registered.
Definition at line 249 of file ad_helpers.cc.
|
protectedinherited |
Register the definition of the index'th dependent variable \(f(\mathbf{X})\).
[in] | index | The index of the entry in the global list of dependent variables that this function belongs to. |
[in] | func | The recorded function that defines a dependent variable. |
Definition at line 678 of file ad_helpers.cc.
|
staticconstexprinherited |
Type definition for the dimension of the associated input and output tensor types.
Definition at line 2647 of file ad_helpers.h.
|
privateinherited |
The independent variables for which we must take into account symmetry when extracting their gradient or Hessian values.
Definition at line 2951 of file ad_helpers.h.
|
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 590 of file ad_helpers.h.
|
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 598 of file ad_helpers.h.
|
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 634 of file ad_helpers.h.
|
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 644 of file ad_helpers.h.
|
protectedinherited |
A list of registered independent variables that have been manipulated for a given set of operations.
Definition at line 650 of file ad_helpers.h.
|
mutableprotectedinherited |
A list of registered independent variables that have been extracted and their sensitivities marked.
Definition at line 656 of file ad_helpers.h.
|
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.
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 748 of file ad_helpers.h.
|
protectedinherited |
A list of registered dependent variables.
Definition at line 753 of file ad_helpers.h.