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::CellLevelBase< 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 | |
EnergyFunctional (const unsigned int n_independent_variables) | |
virtual | ~EnergyFunctional ()=default |
Dependent variables | |
void | register_energy_functional (const ad_type &energy) |
scalar_type | compute_energy () const |
void | compute_residual (Vector< scalar_type > &residual) const override |
virtual void | compute_linearization (FullMatrix< scalar_type > &linearization) const override |
Public Member Functions inherited from Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType > | |
CellLevelBase (const unsigned int n_independent_variables, const unsigned int n_dependent_variables) | |
virtual | ~CellLevelBase ()=default |
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 |
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) |
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 helper class that facilitates the implementation of a generic (incremental) variational formulation from which the computation of the residual vector, as well as its linearization, is automated. This class would typically be used to derive the residual vector and tangent matrix (defined on the level of a cell), or a linearized system of equations, starting from a scalar energy functional.
An example of its usage in the case of a residual and tangent computations might be as follows (in this case we'll compute the linearization of a finite-strain hyperelastic solid from a stored/strain energy density function):
In most use cases, and in particular in the code example shown above, the number of independent variables equals the number of dofs_per_cell
for the used finite element.
Definition at line 1224 of file ad_helpers.h.
using Differentiation::AD::EnergyFunctional< 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 1232 of file ad_helpers.h.
using Differentiation::AD::EnergyFunctional< 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 1239 of file ad_helpers.h.
Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >::EnergyFunctional | ( | 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 \(\Psi(\mathbf{X})\), this will be the number of inputs \(\mathbf{X}\), i.e., the dimension of the domain space. |
Definition at line 802 of file ad_helpers.cc.
|
virtualdefault |
Destructor
void Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >::register_energy_functional | ( | const ad_type & | energy | ) |
Register the definition of the total cell energy \(\Psi(\mathbf{X})\).
[in] | energy | A recorded function that defines the total cell energy. This represents the single dependent variable from which both the residual and its linearization are to be computed. |
Definition at line 811 of file ad_helpers.cc.
EnergyFunctional< ADNumberTypeCode, ScalarType >::scalar_type Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >::compute_energy | ( | ) | const |
Evaluation of the total scalar energy functional for a chosen set of degree of freedom values, i.e.
\[ \Psi(\mathbf{X}) \vert_{\mathbf{X}} \]
The values at the evaluation point \(\mathbf{X}\) are obtained by calling CellLevelBase::set_dof_values().
Definition at line 823 of file ad_helpers.cc.
|
overridevirtual |
Evaluation of the residual for a chosen set of degree of freedom values. Underlying this is the computation of the gradient (first derivative) of the scalar function \(\Psi\) with respect to all independent variables, i.e.
\[ \mathbf{r}(\mathbf{X}) = \frac{\partial\Psi(\mathbf{X})}{\partial\mathbf{X}} \Big\vert_{\mathbf{X}} \]
The values at the evaluation point \(\mathbf{X}\) are obtained by calling CellLevelBase::set_dof_values().
[out] | residual | A Vector object, for which the value for each entry represents the residual value for the corresponding local degree of freedom. The output residual vector has a length corresponding to n_independent_variables . |
Implements Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >.
Definition at line 877 of file ad_helpers.cc.
|
overridevirtual |
Compute the linearization of the residual vector around a chosen set of degree of freedom values. Underlying this is the computation of the Hessian (second derivative) of the scalar function \(\Psi\) with respect to all independent variables, i.e.
\[ \frac{\partial\mathbf{r}(\mathbf{X})}{\partial\mathbf{X}} = \frac{\partial^{2}\Psi(\mathbf{X})}{\partial\mathbf{X} \otimes \partial\mathbf{X}} \Big\vert_{\mathbf{X}} \]
The values at the evaluation point \(\mathbf{X}\) are obtained by calling CellLevelBase::set_dof_values().
[out] | linearization | A FullMatrix representing the linearization of the residual vector. The output linearization matrix has dimensions corresponding to n_independent_variables \(\times\)n_independent_variables . |
Implements Differentiation::AD::CellLevelBase< ADNumberTypeCode, ScalarType >.
Definition at line 942 of file ad_helpers.cc.