Reference documentation for deal.II version 9.0.0
|
#include <deal.II/lac/solver_fire.h>
Classes | |
struct | AdditionalData |
Public Member Functions | |
SolverFIRE (SolverControl &solver_control, VectorMemory< VectorType > &vector_memory, const AdditionalData &data=AdditionalData()) | |
SolverFIRE (SolverControl &solver_control, const AdditionalData &data) | |
virtual | ~SolverFIRE () |
template<typename PreconditionerType = DiagonalMatrix<VectorType>> | |
void | solve (const std::function< double(VectorType &, const VectorType &)> &compute, VectorType &x, const PreconditionerType &inverse_mass_matrix) |
template<typename MatrixType , typename PreconditionerType > | |
void | solve (const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner) |
Public Member Functions inherited from Solver< VectorType > | |
Solver (SolverControl &solver_control, VectorMemory< VectorType > &vector_memory) | |
Solver (SolverControl &solver_control) | |
boost::signals2::connection | connect (const std::function< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType ¤t_iterate)> &slot) |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Protected Member Functions | |
virtual void | print_vectors (const unsigned int, const VectorType &x, const VectorType &v, const VectorType &g) const |
Protected Attributes | |
const AdditionalData | additional_data |
Protected Attributes inherited from Solver< VectorType > | |
GrowingVectorMemory< VectorType > | static_vector_memory |
VectorMemory< VectorType > & | memory |
boost::signals2::signal< SolverControl::State(const unsigned int iteration, const double check_value, const VectorType ¤t_iterate), StateCombiner > | iteration_status |
Additional Inherited Members | |
Public Types inherited from Solver< VectorType > | |
typedef VectorType | vector_type |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
FIRE (Fast Inertial Relaxation Engine) for minimization of (potentially non-linear) objective function \(E(\mathbf x)\), \(\mathbf x\) is a vector of \(n\) variables ( \(n\) is the number of variables of the objective function). Like all other solver classes, it can work on any kind of vector and matrix as long as they satisfy certain requirements (for the requirements on matrices and vectors in order to work with this class, see the documentation of the Solver base class). The type of the solution vector must be passed as template argument, and defaults to Vector<double>.
FIRE is a damped dynamics method described in Structural Relaxation Made Simple by Bitzek et al. 2006, typically used to find stable equilibrium configurations of atomistic systems in computational material science. Starting from a given initial configuration of the atomistic system, the algorithm relies on inertia to obtain (nearest) configuration with least potential energy.
Notation:
Given initial values for \(\Delta t\), \(\alpha = \alpha_0\), \(\epsilon\), \(\mathbf x = \mathbf x_0\) and \(\mathbf v= \mathbf 0\) along with a given mass matrix \(\mathbf M\), FIRE algorithm is as follows,
Also see Energy-Minimization in Atomic-to-Continuum Scale-Bridging Methods by Eidel et al. 2011.
Definition at line 89 of file solver_fire.h.
SolverFIRE< VectorType >::SolverFIRE | ( | SolverControl & | solver_control, |
VectorMemory< VectorType > & | vector_memory, | ||
const AdditionalData & | data = AdditionalData() |
||
) |
Constructor.
SolverFIRE< VectorType >::SolverFIRE | ( | SolverControl & | solver_control, |
const AdditionalData & | data | ||
) |
Constructor. Use an object of type GrowingVectorMemory as a default to allocate memory.
|
virtual |
Virtual destructor.
void SolverFIRE< VectorType >::solve | ( | const std::function< double(VectorType &, const VectorType &)> & | compute, |
VectorType & | x, | ||
const PreconditionerType & | inverse_mass_matrix | ||
) |
Obtain a set of variables x
that minimize an objective function described by the polymorphic function wrapper compute
, with a given preconditioner inverse_mass_matrix
and initial x
values. The function compute
returns the objective function's value and updates the objective function's gradient (with respect to the variables) when passed in as first argument based on the second argument– the state of variables.
void SolverFIRE< VectorType >::solve | ( | const MatrixType & | A, |
VectorType & | x, | ||
const VectorType & | b, | ||
const PreconditionerType & | preconditioner | ||
) |
Solve for x that minimizes \(E(\mathbf x)\) for the special case when \(E(\mathbf x) = \frac{1}{2} \mathbf x^{T} \mathbf A \mathbf x - \mathbf x^{T} \mathbf b\).
|
protectedvirtual |
Interface for derived class. This function gets the current iteration x
(variables), v
(x's time derivative) and g
(the gradient) in each step. It can be used for graphical output of the convergence history.
|
protected |
Additional data to the solver.
Definition at line 187 of file solver_fire.h.