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::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 | |
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 |
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 |
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 () |
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).
Definition at line 840 of file ad_helpers.h.
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 848 of file ad_helpers.h.
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 855 of file ad_helpers.h.
Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >::CellLevelBase | ( | 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 719 of file ad_helpers.cc.
|
virtualdefault |
Destructor
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.
[in] | dof_values | A 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. |
scalar_type
. Depending on the selected ADNumberTypeCode
, this may or may not correspond with the ScalarType
prescribed as a template argument.Definition at line 730 of file ad_helpers.cc.
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.
[in] | values | A 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_indices | A 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() . |
const std::vector< typename CellLevelBase< ADNumberTypeCode, ScalarType >::ad_type > & Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >::get_sensitive_dof_values | ( | ) | const |
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().
Definition at line 753 of file ad_helpers.cc.
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.
[in] | dof_values | A vector field associated with local degree of freedom values on the current finite element. These define 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_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 778 of file ad_helpers.cc.
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.
[in] | values | A vector field from which the values of all independent variables is to be extracted. |
[in] | local_dof_indices | A 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() . |
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.
|
pure virtual |
Compute the value of the residual vector field \(\mathbf{r}(\mathbf{X})\).
[out] | residual | A Vector object with the value for each component of the vector field evaluated at the point defined by the independent variable values. |
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::ResidualLinearization< ADNumberTypeCode, ScalarType >, and Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >.
|
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}} \]
[out] | linearization | A FullMatrix with the gradient of each component of the vector field evaluated at the point defined by the independent variable values. |
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::ResidualLinearization< ADNumberTypeCode, ScalarType >, and Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >.