Reference documentation for deal.II version 9.1.1
|
#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 Types inherited from Differentiation::AD::PointLevelFunctionsBase< dim, ADNumberTypeCode, ScalarType > | |
using | scalar_type = typename HelperBase< ADNumberTypeCode, ScalarType >::scalar_type |
using | ad_type = typename HelperBase< ADNumberTypeCode, ScalarType >::ad_type |
Public Types inherited from Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType > | |
using | scalar_type = typename AD::NumberTraits< ScalarType, ADNumberTypeCode >::scalar_type |
using | ad_type = typename AD::NumberTraits< ScalarType, ADNumberTypeCode >::ad_type |
Public Member Functions | |
Constructor / destructor | |
ScalarFunction (const unsigned int n_independent_variables) | |
virtual | ~ScalarFunction ()=default |
Public Member Functions inherited from Differentiation::AD::PointLevelFunctionsBase< dim, ADNumberTypeCode, ScalarType > | |
PointLevelFunctionsBase (const unsigned int n_independent_variables, const unsigned int n_dependent_variables) | |
virtual | ~PointLevelFunctionsBase ()=default |
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) |
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< ad_type > | get_sensitive_variables (const ExtractorType &extractor) const |
Public Member Functions inherited from Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType > | |
HelperBase (const unsigned int n_independent_variables, const unsigned int n_dependent_variables) | |
virtual | ~HelperBase ()=default |
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 |
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 () |
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) |
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 3094 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 3103 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 3110 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 1330 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 1344 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 1357 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 1407 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 1481 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 1580 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 1617 of file ad_helpers.cc.