#include <energyMinimizer.h>
Inheritance diagram for EnergyMinimizer:
Public Member Functions | |
bool | operator== (const EnergyMinimizer &energy_minimizer) throw () |
Equality operator. | |
Constructors and Destructors | |
EnergyMinimizer () | |
Default constructor. | |
EnergyMinimizer (ForceField &force_field) | |
Constructor. | |
EnergyMinimizer (ForceField &force_field, const Options &options) | |
Constructor. | |
EnergyMinimizer (const EnergyMinimizer &energy_minimizer) | |
Copy constructor. | |
virtual | ~EnergyMinimizer () |
Destructor. | |
Assignments | |
const EnergyMinimizer & | operator= (const EnergyMinimizer &energy_minimizer) |
Assignment operator. | |
Debugging and Diagnostics | |
bool | isValid () const |
Is the energy minimizer valid : did the setup work? | |
Setup methods | |
bool | setup (ForceField &force_field) |
Sets up the energy minimizer. | |
bool | setup (ForceField &force_field, SnapShotManager *ssm) |
Sets up the energy minimizer. | |
bool | setup (ForceField &force_field, SnapShotManager *ssm, const Options &options) |
Sets up the energy minimizer. | |
bool | setup (ForceField &force_field, const Options &options) |
Sets up the energy minimizer. | |
virtual bool | specificSetup () |
Specific setup. | |
Accessors | |
virtual bool | isConverged () const |
Implements the convergence criterion. | |
virtual double | findStep () |
Calculate the next step. | |
virtual void | updateDirection () |
Update the search direction. | |
virtual double | updateEnergy () |
Update energy. | |
virtual void | updateForces () |
Update forces and store them in current_grad_. | |
void | storeGradientEnergy () |
Store the current energy and gradient. | |
virtual void | printEnergy () const |
Print the energy. | |
virtual void | takeSnapShot () const |
Take a snapshot of the system. | |
virtual void | finishIteration () |
Finishing step for this iteration. | |
Size | getNumberOfIterations () const |
Return the number of iterations performed. | |
Gradient & | getDirection () |
Return a reference to the current search direction. | |
Gradient & | getGradient () |
Return a reference to the current gradient. | |
Gradient & | getInitialGradient () |
Return a reference to the initial gradient. | |
double | getEnergy () const |
Return the current energy. | |
double & | getEnergy () |
Return a reference to the current energy. | |
double | getInitialEnergy () const |
Return the initial energy. | |
double & | getInitialEnergy () |
Return a mutable reference to the initial energy. | |
void | setNumberOfIterations (Size number_of_iterations) |
Set the number of iterations performed so far. | |
Size | getMaxNumberOfIterations () const |
Get the maximum number of iterations. | |
void | setMaxNumberOfIterations (Size number_of_iterations) |
Set the maximum number of iterations. | |
void | setMaxSameEnergy (Size number) |
Set the maximum number of iterations allowed with equal energy (second convergence criterion). | |
Size | getMaxSameEnergy () const |
Get the maximum number of iterations allowed with equal energy (second convergence criterion). | |
void | setEnergyOutputFrequency (Size energy_output_frequency) |
Set the energy output frequency. | |
Size | getEnergyOutputFrequency () const |
Get the energy ouput frequency. | |
void | setEnergyDifferenceBound (float energy_difference_bound) |
Set the energy difference bound for convergence. | |
float | getEnergyDifferenceBound () const |
Get the energy difference bound. | |
void | setMaxGradient (float max_gradient) |
Set the maximum RMS gradient (first convergence criterion). | |
float | getMaxGradient () const |
Get the maximum RMS gradient (first convergence criterion). | |
void | setMaximumDisplacement (float maximum_displacement) |
Set the maximum displacement value. | |
float | getMaximumDisplacement () const |
Get the maximum displacement value. | |
void | setSnapShotFrequency (Size snapshot_frequency) |
Set the snapshot frequency. | |
Size | getSnapShotFrequency () const |
Get the snapshot output frequency. | |
ForceField * | getForceField () |
Return the force field of the energy minimizer. | |
Size | getForceUpdateCounter () const throw () |
Return the number of force updates since the start of the minimization. | |
Size | getEnergyUpdateCounter () const throw () |
Return the number of energy updates since the start of the minimization. | |
virtual bool | minimize (Size steps=0, bool resume=false) |
Minimize the energy of the system bound to the force field. | |
void | enableEnergyAbortCondition (bool state) |
Specify if the MDSimulation aborts if the Energy is greater than abort_energy_. | |
bool | energyAbortConditionEnabled () const |
Query if the MDSimulation aborts if the Energy is greater than abort_energy_. | |
void | setEnergyToAbort (float value) |
Set the value for the energy, that will result in aborting the minization, if it will be surpassed. | |
float | getEnergyToAbort () const |
bool | wasAborted () const throw () |
Return true, if the minimization was aborted, e.g. | |
Public Attributes | |
Public Attributes | |
Options | options |
Options. | |
Protected Attributes | |
Protected Attributes | |
Gradient | initial_grad_ |
The gradient at the beginning of the current minimization step. | |
Gradient | current_grad_ |
The current gradient. | |
double | initial_energy_ |
The energy at the beginning of the current minimization step. | |
double | current_energy_ |
The current energy. | |
Gradient | old_grad_ |
The gradient from the last step. | |
double | old_energy_ |
The energy from the last step. | |
Gradient | direction_ |
The current search direction. | |
bool | valid_ |
The boolean variable indicates if the setup of the energy minimizer was successful. | |
SnapShotManager * | snapshot_ |
Pointer to a SnapShotManager for storing snapshots of the system. | |
ForceField * | force_field_ |
The force field bound to the energy minimizer. | |
Size | number_of_iterations_ |
The current iteration number. | |
Size | maximal_number_of_iterations_ |
Maximum number of iterations. | |
Size | energy_output_frequency_ |
Frequency of energy output. | |
Size | snapshot_frequency_ |
Frequency of atom coordinate ouput. | |
double | energy_difference_bound_ |
If the energy difference (before and after an iteration) is smaller than this bound, the minimization procedure stops. | |
double | max_gradient_ |
The maximum RMS gradient tolerated (first convergence criterion). | |
Size | max_same_energy_ |
The maximum number of iterations with same energy. | |
Size | same_energy_counter_ |
A counter for the number of steps with a similar energy. | |
float | maximum_displacement_ |
The maximal shift of an atom per iteration step (in Angstrom). | |
Size | force_update_counter_ |
Internal counter: how often is a force update done. | |
Size | energy_update_counter_ |
Internal counter: how often is an energy update done. | |
float | cutlo_ |
Numerical lower bound: we don't want to compute the reciprocal of a number which is lower than 'cutlo_'. | |
double | step_ |
The last step size (in respect of the length of the computed direction vector), so the length of the last step was . | |
bool | abort_by_energy_enabled_ |
float | abort_energy_ |
bool | aborted_ |
Base class for all minimizer for geometry optimization.
|
Default constructor.
|
|
Constructor.
|
|
Constructor.
|
|
Destructor.
|
|
Calculate the next step. This method is implemented in each minimizer class and tries to determine the next step to be taken. It typically performs a line search. The value returned is usually the step length with respect to the current direction.
Reimplemented in ConjugateGradientMinimizer, ShiftedLVMMMinimizer, SteepestDescentMinimizer, and StrangLBFGSMinimizer. |
|
Finishing step for this iteration. This method should be called at the end of the main iteration loop implemented in minimize . It takes over some administrative stuff:
This method should be overwritten only in rare cases. Even then, the programmer should make sure to call All derived classes should call this method at the end of the minimize main loop. Otherwise strange things might happen.
|
|
Return the number of energy updates since the start of the minimization.
|
|
Return the number of force updates since the start of the minimization.
|
|
Get the maximum RMS gradient (first convergence criterion). The gradient unit of the gradient is kJ/(mol ). |
|
Return the number of iterations performed.
|
|
Implements the convergence criterion. If the convergence criterion is fulfilled, this method returns true. The convergence criterion is implemented as one of two conditions: (1) {RMS gradient} is below max_rms_gradient_ (2) same_energy_counter_ is above max_same_energy_ If any of these conditions hold isConverged returns true. This method should be reimplemented in derived classes for a different convergence criterion. |
|
Minimize the energy of the system bound to the force field.
If a number of steps is given, the minimization is aborted after that number of steps, regardless of the number of steps given in the options (
Reimplemented in ConjugateGradientMinimizer, ShiftedLVMMMinimizer, SteepestDescentMinimizer, and StrangLBFGSMinimizer. |
|
Print the energy.
This method is called by finishIteration after every energy_output_frequency_ steps. It prints the current RMS gradient and the current energy to Log |
|
Set the value for the energy, that will result in aborting the minization, if it will be surpassed. Default value: 10^9. |
|
Set the maximum RMS gradient (first convergence criterion). The gradient unit of the gradient is kJ/(mol ). |
|
Set the maximum displacement value. This is the maximum distance an atom may be moved by the minimizer in one iteration. |
|
Set the number of iterations performed so far.
|
|
Sets up the energy minimizer.
|
|
Sets up the energy minimizer.
|
|
Sets up the energy minimizer.
|
|
Sets up the energy minimizer.
|
|
Store the current energy and gradient. The current gradient and current energy is copied into initial energy and initial gradient. This is usually done at the start of an iteration. |
|
Take a snapshot of the system.
This method is called by finishIteration after every snapshot_frequency_ steps. It saves a
|
|
Update the search direction. This method is implemented by the derived classes to implement a method to determine a new search direction. Reimplemented in ConjugateGradientMinimizer, ShiftedLVMMMinimizer, SteepestDescentMinimizer, and StrangLBFGSMinimizer. |
|
Update energy.
This method calls |
|
Update forces and store them in current_grad_.
This method calls |
|
Return true, if the minimization was aborted, e.g. because of strange energies or gradient. |
|
The current energy.
|
|
The current gradient.
|
|
Numerical lower bound: we don't want to compute the reciprocal of a number which is lower than 'cutlo_'.
|
|
If the energy difference (before and after an iteration) is smaller than this bound, the minimization procedure stops.
|
|
Internal counter: how often is an energy update done. Measure for the speed of minimization. |
|
The force field bound to the energy minimizer. Among other data the force field contains the molecular system whose energy will be minimized by the energy minimizer. |
|
Internal counter: how often is a force update done. Measure for the speed of minimization. |
|
The energy at the beginning of the current minimization step.
|
|
The gradient at the beginning of the current minimization step.
|
|
The maximum number of iterations with same energy. When this number is reached, we assume the system to have converged (second convergence criterion) |
|
The maximal shift of an atom per iteration step (in Angstrom).
|
|
A counter for the number of steps with a similar energy.
|
|
The last step size (in respect of the length of the computed direction vector), so the length of the last step was .
|