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 | |
VectorFunction (const unsigned int n_independent_variables, const unsigned int n_dependent_variables) | |
virtual | ~VectorFunction ()=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_variables (const std::vector< ad_type > &funcs) |
template<typename ValueType , typename ExtractorType > | |
void | register_dependent_variable (const ValueType &funcs, const ExtractorType &extractor) |
void | compute_values (Vector< scalar_type > &values) const |
void | compute_jacobian (FullMatrix< scalar_type > &jacobian) const |
template<typename ExtractorType_Row > | |
static internal::VectorFieldValue< dim, scalar_type, ExtractorType_Row >::type | extract_value_component (const Vector< scalar_type > &values, const ExtractorType_Row &extractor_row) |
template<typename ExtractorType_Row , typename ExtractorType_Col > | |
static internal::VectorFieldJacobian< dim, scalar_type, ExtractorType_Row, ExtractorType_Col >::type | extract_jacobian_component (const FullMatrix< scalar_type > &jacobian, const ExtractorType_Row &extractor_row, const ExtractorType_Col &extractor_col) |
static Tensor< 0, dim, scalar_type > | extract_jacobian_component (const FullMatrix< scalar_type > &jacobian, const FEValuesExtractors::Scalar &extractor_row, const FEValuesExtractors::Scalar &extractor_col) |
static SymmetricTensor< 4, dim, scalar_type > | extract_jacobian_component (const FullMatrix< scalar_type > &jacobian, const FEValuesExtractors::SymmetricTensor< 2 > &extractor_row, const FEValuesExtractors::SymmetricTensor< 2 > &extractor_col) |
A helper class that facilitates the evaluation of a vector of functions, typically one that represents a collection of coupled, multi-dimensional fields. This class would typically be used to compute the linearization of a set of kinetic field variables defined at the quadrature point level.
An example of its usage in the case of linearizing the kinetic variables derived from a multi-field constitutive law might be as follows:
Definition at line 3487 of file ad_helpers.h.
using Differentiation::AD::VectorFunction< 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 3496 of file ad_helpers.h.
using Differentiation::AD::VectorFunction< 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 3503 of file ad_helpers.h.
Differentiation::AD::VectorFunction< dim, ADNumberTypeCode, ScalarType >::VectorFunction | ( | const unsigned int | n_independent_variables, |
const unsigned int | n_dependent_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. |
[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. |
Definition at line 1655 of file ad_helpers.cc.
|
virtualdefault |
Destructor.
void Differentiation::AD::VectorFunction< dim, ADNumberTypeCode, ScalarType >::register_dependent_variables | ( | const std::vector< ad_type > & | funcs | ) |
Register the definition of the vector field \(\boldsymbol{\Psi}(\mathbf{X})\).
[in] | funcs | A vector of recorded functions that defines the dependent variables. |
Definition at line 1670 of file ad_helpers.cc.
void Differentiation::AD::VectorFunction< dim, ADNumberTypeCode, ScalarType >::register_dependent_variable | ( | const ValueType & | funcs, |
const ExtractorType & | extractor | ||
) |
Register the definition of the vector field \(\hat{\mathbf{g}}(\mathbf{X}) \subset \boldsymbol{\Psi}(\mathbf{X})\) that may represent a subset of the dependent variables.
[in] | funcs | The recorded functions that define a set of dependent variables. |
[in] | extractor | An extractor associated with the input field variables. This effectively defines which components of the global set of dependent variables this field is associated with. |
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>.void Differentiation::AD::VectorFunction< dim, ADNumberTypeCode, ScalarType >::compute_values | ( | Vector< scalar_type > & | values | ) | const |
Compute the value of the vector field \(\boldsymbol{\Psi}(\mathbf{X})\).
[out] | values | A Vector object with the value for each component of the vector field evaluated at the point defined by the independent variable values. The output values vector has a length corresponding to n_dependent_variables . |
Definition at line 1686 of file ad_helpers.cc.
void Differentiation::AD::VectorFunction< dim, ADNumberTypeCode, ScalarType >::compute_jacobian | ( | FullMatrix< scalar_type > & | jacobian | ) | const |
Compute the Jacobian (first derivative) of the vector field with respect to all independent variables, i.e.
\[ \mathbf{J}(\boldsymbol{\Psi}) = \frac{\partial\boldsymbol{\Psi}(\mathbf{X})}{\partial\mathbf{X}} \]
[out] | jacobian | A FullMatrix with the gradient of each component of the vector field evaluated at the point defined by the independent variable values. The output jacobian matrix has dimensions corresponding to n_dependent_variables \(\times\)n_independent_variables . |
Definition at line 1742 of file ad_helpers.cc.
|
static |
Extract the set of functions' values for a subset of dependent variables \(\mathbf{g} \subset \boldsymbol{\Psi}(\mathbf{X})\).
[in] | values | A Vector object with the value for each component of the vector field evaluated at the point defined by the independent variable values. |
[in] | extractor_row | An extractor associated with the input field variables. This effectively defines which components of the global set of dependent variables this field is associated with. |
|
static |
Extract the Jacobian of the subset of dependent functions \(\mathbf{g} \subset \boldsymbol{\Psi}(\mathbf{X})\) for a subset of independent variables \(\mathbf{A} \subset \mathbf{X}\), i.e.
\[ \mathbf{J}(\mathbf{g}) = \frac{\partial\mathbf{g}(\mathbf{X})}{\partial\mathbf{A}} \]
The first index of the Jacobian matrix \(\mathbf{J}(\mathbf{g})\) relates to the dependent variables, while the second index relates to the independent variables.
[in] | jacobian | The Jacobian of the vector function with respect to all independent variables, i.e., that returned by compute_jacobian(). |
[in] | extractor_row | An extractor associated with the input field variables for which the first index of the Jacobian is extracted. This effectively defines the correspondence between components of the global set of dependent variables and the field (representing a subset of dependent functions) associated with the extractor. |
[in] | extractor_col | An extractor associated with the input field variables for which the second index of the Jacobian is extracted. This effectively defines the correspondence between components of the global set of independent variables and the field (representing a subset of independent variables) associated with the extractor. |
extractor_row
and extractor_col
. This corresponds to subsetting a whole set of rows and columns of the Jacobian 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 Jacobian of the subset of dependent functions \(\mathbf{g} \subset \boldsymbol{\Psi}(\mathbf{X})\) for a subset of independent variables \(\mathbf{A} \subset \mathbf{X}\), i.e.
\[ \mathbf{J}(\mathbf{g}) = \frac{\partial\mathbf{g}(\mathbf{X})}{\partial\mathbf{A}} \]
This function is a specialization of the above for rank-0 tensors (scalars). This corresponds to extracting a single entry of the Jacobian matrix because both extractors imply selection of just a single row or column of the matrix.
Definition at line 1821 of file ad_helpers.cc.
|
static |
Extract the Jacobian of the subset of dependent functions \(\mathbf{g} \subset \boldsymbol{\Psi}(\mathbf{X})\) for a subset of independent variables \(\mathbf{A} \subset \mathbf{X}\), i.e.
\[ \mathbf{J}(\mathbf{g}) = \frac{\partial\mathbf{g}(\mathbf{X})}{\partial\mathbf{A}} \]
This function is a specialization of the above for rank-4 symmetric tensors.
Definition at line 1859 of file ad_helpers.cc.