Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Types | List of all members
Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType > Class Template Reference

#include <deal.II/differentiation/ad/ad_helpers.h>

Inheritance diagram for Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >:
[legend]

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 ()
 

Additional Inherited Members

- Static Public Member Functions inherited from Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >
static void configure_tapeless_mode (const unsigned int n_independent_variables, const bool ensure_persistent_setting=true)
 
- Protected Member Functions inherited from Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >
void activate_tape (const typename Types< ad_type >::tape_index tape_index, const bool read_mode)
 
void reset_registered_independent_variables ()
 
void set_sensitivity_value (const unsigned int index, const scalar_type &value)
 
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
 
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)
 
- Protected Attributes inherited from Differentiation::AD::HelperBase< ADNumberTypeCode, ScalarType >
TapedDrivers< ad_type, scalar_typetaped_driver
 
TapelessDrivers< ad_type, scalar_typetapeless_driver
 
std::vector< scalar_typeindependent_variable_values
 
std::vector< ad_typeindependent_variables
 
std::vector< bool > registered_independent_variable_values
 
std::vector< bool > registered_marked_independent_variables
 
std::vector< ad_typedependent_variables
 
std::vector< bool > registered_marked_dependent_variables
 

Detailed Description

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
class Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >

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):

// Existing data structures:
Vector<double> solution (...); // Or another vector type
std::vector<types::global_dof_index> local_dof_indices (...);
const FEValuesExtractors::Vector u_fe (...);
FEValues<dim> fe_values (...);
const unsigned int n_q_points (...);
Vector<double> cell_rhs (...);
// Assembly loop:
for (auto & cell : ...)
{
cell->get_dof_indices(local_dof_indices);
const unsigned int n_independent_variables =
local_dof_indices.size();
// Create some aliases for the AD helper.
// In the example, the AD_typecode used for the template argument can
// be refer to either a taped or tapeless type.
using ADHelper = AD::EnergyFunctional<...>;
using ADNumberType = typename ADHelper::ad_type;
// Create and initialize an instance of the helper class.
ADHelper ad_helper(n_independent_variables);
// Initialize the local data structures for assembly.
// This is also taken care of by the ADHelper, so this step could
// be skipped.
cell_rhs.reinit(n_independent_variables);
// An optional call to set the amount of memory to be allocated to
// storing taped data.
// If using a taped AD number then we would likely want to increase
// the buffer size from the default values as the expression for each
// residual component will likely be lengthy, and therefore memory
// intensive.
ad_helper.set_tape_buffer_sizes(...);
// If using a taped AD number, then at this point we would initiate
// taping of the expression for the energy for this FE type and
// material combination:
// Select a tape number to record to
const typename Types<ad_type>::tape_index tape_index = ...;
// Indicate that we are about to start tracing the operations for
// function evaluation on the tape. If this tape has already been
// used (i.e. the operations are already recorded) then we
// (optionally) load the tape and reuse this data.
const bool is_recording
= ad_helper.start_recording_operations(tape_index);
// The steps that follow in the recording phase are required for
// tapeless methods as well.
if (is_recording == true)
{
// This is the "recording" phase of the operations.
// First, we set the values for all DoFs.
ad_helper.register_dof_values(solution, local_dof_indices);
// Then we get the complete set of degree of freedom values as
// represented by auto-differentiable numbers. The operations
// performed with these variables are tracked by the AD library
// from this point until stop_recording_operations() is called.
const std::vector<ADNumberType> dof_values_ad
= ad_helper.get_sensitive_dof_values();
// Then we do some problem specific tasks, the first being to
// compute all values, gradients, etc. based on sensitive AD DoF
// values. Here we are fetching the displacement gradients at each
// quadrature point.
std::vector<Tensor<2, dim, ADNumberType>> Grad_u(
fe_values[u_fe].get_function_gradients_from_local_dof_values(
dof_values_ad, Grad_u);
// This variable stores the cell total energy.
// IMPORTANT: Note that it is hand-initialized with a value of
// zero. This is a highly recommended practise, as some AD numbers
// appear not to safely initialize their internal data structures.
ADNumberType energy_ad = ADNumberType(0.0);
// Compute the cell total energy = (internal + external) energies
for (unsigned int q_point = 0; q_point < n_q_points; ++q_point)
{
// Calculate the deformation gradient at this quadrature point
unit_symmetric_tensor<dim>() + Grad_u[q_point];
ExcMessage("Negative determinant of the deformation "
"gradient detected!"));
// Add contribution of the internal energy:
// Integrate the stored energy density function with the current
// solution.
energy_ad += get_Psi(F) * fe_values.JxW(q_point);
}
// Add contribution from external energy:
// Loop over faces and accumulate external energy into cell
// total energy.
for (unsigned int face : ...)
if (cell->face(face)->at_boundary())
energy_ad += ...
// Register the definition of the total cell energy
ad_helper.register_energy_functional(energy_ad);
// Indicate that we have completed tracing the operations onto
// the tape.
ad_helper.stop_recording_operations(false); // write_tapes_to_file
}
else
{
// This is the "tape reuse" phase of the operations.
// Here we will leverage the already traced operations that reside
// on a tape, and simply re-evaluate the tape at a different point
// to get the function values and their derivatives.
// Load the existing tape to be reused
ad_helper.activate_recorded_tape(tape_index);
// Set the new values of the independent variables where the
// recorded dependent functions are to be evaluated (and
// differentiated around).
ad_helper.set_dof_values(solution, local_dof_indices);
}
// Compute the residual values and their Jacobian at the
// evaluation point
ad_helper.compute_residual(cell_rhs);
cell_rhs *= -1.0; // RHS = - residual
ad_helper.compute_linearization(cell_matrix);
}

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.

Warning
ADOL-C does not support the standard threading models used by deal.II, so this class should not be embedded within a multithreaded function when using ADOL-C number types. It is, however, suitable for use in both serial and MPI routines.
Author
Jean-Paul Pelteret, 2016, 2017, 2018

Definition at line 1224 of file ad_helpers.h.

Member Typedef Documentation

◆ scalar_type

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
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.

◆ ad_type

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
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.

Constructor & Destructor Documentation

◆ EnergyFunctional()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >::EnergyFunctional ( const unsigned int  n_independent_variables)

The constructor for the class.

Parameters
[in]n_independent_variablesThe 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.
Note
There is only one dependent variable associated with the total energy attributed to the local finite element. That is to say, this class assumes that the (local) right hand side and matrix contribution is computed from the first and second derivatives of a scalar function \(\Psi(\mathbf{X})\).

Definition at line 802 of file ad_helpers.cc.

◆ ~EnergyFunctional()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType = double>
virtual Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >::~EnergyFunctional ( )
virtualdefault

Destructor

Member Function Documentation

◆ register_energy_functional()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >::register_energy_functional ( const ad_type energy)

Register the definition of the total cell energy \(\Psi(\mathbf{X})\).

Parameters
[in]energyA 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.
Note
For this class that expects only a single scalar dependent variable, this function must only be called once per tape.
For taped AD numbers, this operation is only valid in recording mode.

Definition at line 811 of file ad_helpers.cc.

◆ compute_energy()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
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().

Returns
The value of the energy functional at the evaluation point corresponding to a chosen set of local degree of freedom values.

Definition at line 823 of file ad_helpers.cc.

◆ compute_residual()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >::compute_residual ( Vector< scalar_type > &  residual) const
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().

Parameters
[out]residualA 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.

◆ compute_linearization()

template<enum AD::NumberTypes ADNumberTypeCode, typename ScalarType >
void Differentiation::AD::EnergyFunctional< ADNumberTypeCode, ScalarType >::compute_linearization ( FullMatrix< scalar_type > &  linearization) const
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().

Parameters
[out]linearizationA 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.


The documentation for this class was generated from the following files: