File Dimer.h

Copyright

This code is licensed under the 3-clause BSD license.

Copyright ETH Zurich, Laboratory for Physical Chemistry, Reiher Group.

See LICENSE.txt for details.

namespace Scine

This header file contains functions that allow for common notation for common things that can be done at a different degree of derivatives.

This header contains alias definitions defining which classes to use for the different degrees of derivatives.

namespace Utils
class Dimer : public Scine::Utils::Optimizer
#include <Dimer.h>

An implementation of the Dimer optimization algorithm for saddlepoints.

transition state search without hessian by finding minimum curvature by gradient calculations of a dimer on PES dimer is constructed by estimated eigenvector or random vector, then dimer is rotated iteratively until it is aligned with lowest curvature mode, then translation is performed this is repeated until converged

Implemented, as described in: Kaestner J. and Sherwood P., Journal of Chemical Physics, 2008 Shang C. and Liu ZP., Journal of Chemical Theory and Computation, 2010

Public Functions

Dimer()

Default constructor.

template<class UpdateFunction>
int optimize(Eigen::VectorXd &parameters, UpdateFunction &&function, GradientBasedCheck &check)

The main routine of the optimizer that carries out the actual optimization.

Return

int Returns the number of optimization cycles carried out until the conclusion of the optimization function.

Template Parameters
  • UpdateFunction: A lambda function with a void return value, and the arguments:

    1. const Eigen::VectorXd& parameters

    2. double& value

    3. Eigen::VectorXd& gradients

Parameters
  • parameters: The parameters to be optimized.

  • function: The function to be evaluated in order to get values and gradients for a given set of parameters.

  • check: The ConvergenceCheck to be used in order to determine when the optimization is finished or should stop for other reasons.

virtual void addSettingsDescriptors(UniversalSettings::DescriptorCollection &collection) const

Adds all relevant options to the given UniversalSettings::DescriptorCollection thus expanding it to include the Dimers’s options.

Parameters
  • collection: The DescriptorCollection to which new fields shall be added.

virtual void applySettings(const Settings &settings)

Updates the Dimer’s options with those values given in the Settings.

Parameters
  • settings: The settings to update the option of the steepest descent with.

Public Members

bool skipFirstRotation = false
bool decreaseRotationGradientThreshold = false
bool gradientInterpolation = false
bool rotationCG = false
bool rotationLBFGS = true
bool onlyOneRotation = false
bool useGdiis = false
bool useProjectionTrustRadius = true
bool multiScale = true
double radius = 0.01
double phiTolerance = 1e-3
double rotationGradientThresholdFirstCycle = 1e-7
double rotationGradientThresholdOtherCycles = 1e-4
double loweredRotationGradientThreshold = 1e-3
double rotationGradientThreshold = 0.0
double gradRMSDthreshold = 1e-3
double trustRadius = 0.5

The maximum RMS of a taken step.

double defaultTranslationStep = 1

The factor to multiple the stepsize vector in a steepest descent step.

unsigned int maxRotationsFirstCycle = 100
unsigned int maxRotationsOtherCycles = 100
unsigned int intervalOfRotations = 5
unsigned int cycleOfRotationGradientDecrease = 5
unsigned int maxBacktracking = 5
std::unique_ptr<Eigen::VectorXd> guessVector

Public Static Attributes

constexpr const char *dimerSkipFirstRotation = "dimer_skip_first_rotation"
constexpr const char *dimerDecreaseRotationGradientThreshold = "dimer_decrease_rotation_gradient_threshold"
constexpr const char *dimerGradientInterpolation = "dimer_gradient_interpolation"
constexpr const char *dimerRotationCG = "dimer_rotation_conjugate_gradient"
constexpr const char *dimerRotationLBFGS = "dimer_rotation_lbfgs"
constexpr const char *dimerOnlyOneRotation = "dimer_only_one_rotation"
constexpr const char *dimerUseGdiis = "dimer_gdiis"
constexpr const char *dimerUseProjectionTrustRadius = "dimer_projection_trust_radius"
constexpr const char *dimerMultiScale = "dimer_multi_scale"
constexpr const char *dimerRadius = "dimer_radius"
constexpr const char *dimerPhiTolerance = "dimer_phi_tolerance"
constexpr const char *dimerRotationGradientThresholdFirstCycle = "dimer_rotation_gradient_first"
constexpr const char *dimerRotationGradientThresholdOtherCycles = "dimer_rotation_gradient_other"
constexpr const char *dimerLoweredRotationGradientThreshold = "dimer_lowered_rotation_gradient"
constexpr const char *dimerGradRMSDthreshold = "dimer_grad_rmsd_threshold"
constexpr const char *dimerTrustRadius = "dimer_trust_radius"
constexpr const char *dimerDefaultTranslationStep = "dimer_default_translation_step"
constexpr const char *dimerMaxRotationsFirstCycle = "dimer_max_rotations_first_cycle"
constexpr const char *dimerMaxRotationsOtherCycles = "dimer_max_rotations_other_cycle"
constexpr const char *dimerIntervalOfRotations = "dimer_interval_of_rotations"
constexpr const char *dimerCycleOfRotationGradientDecrease = "dimer_cycle_of_rotation_gradient_decrease"
constexpr const char *dimerMaxBacktracking = "dimer_max_backtracking"

Private Functions

void setVectorsToZero(const int &nParams)
void sanityCheck()
void createDimerAxis(const int &nParams)

Create first dimer.

If guessvector is given, it is used for the dimer axis if no guess, a random vector is used

Return

void

Parameters
  • nParams: The number of parameters to be optimized

template<class UpdateFunction>
bool determineIfPerformRotation(const int &cycle, const Eigen::VectorXd &parameters, UpdateFunction &&function)

Determine whether a rotation of the dimer shall be carried out in this optimization step.

Return

bool Returns whether rotation shall be carried out

Parameters
  • cycle: The number of the current cycle

  • parameters: The parameters to be optimized

Template Parameters
  • UpdateFunction: A lambda function with a void return value, and the arguments:

    1. const Eigen::VectorXd& parameters

    2. double& value

    3. Eigen::VectorXd& gradients

Eigen::VectorXd conjugateGradientRotation(const int &rotationCycles)

Determine step vector of rotation with the Conjugate Gradient (Polak Ribiere) method.

Return

Eigen::VectorXd Returns step vector of rotation

Parameters
  • rotationCycles: The number of individual rotation in the current rotation - 1

Eigen::VectorXd lbfgsRotation(const int &rotationCycles, unsigned int &m)

Determine the step vector of rotation based on the L-BFGS methods.

Return

Eigen::VectorXd Returns step vector of rotation

Parameters
  • rotationCycles: The number of individual rotation in the current rotation - 1

  • m: The counter for the L-BFGS backtracking

template<class UpdateFunction>
void determineDirectionOfDimerToBeMaximized(const double &value, const Eigen::VectorXd &parameters, double &value1, double &curvature, UpdateFunction &&function)

Flips dimeraxis if energy decreases along axis.

Rotate dimer axis if necessary, so dimer axis is afterwards correctly aligned with the direction in which the mode can be optimized

Return

void

Parameters
  • value: The value of the optimized point

  • parameters: The parameters to be optimized

  • value1: The value at R1

  • curvature: The curvature along the dimer axis

Template Parameters
  • UpdateFunction: A lambda function with a void return value, and the arguments:

    1. const Eigen::VectorXd& parameters

    2. double& value

    3. Eigen::VectorXd& gradients

template<class UpdateFunction>
double rotationWithPhi(const Eigen::VectorXd &parameters, const double &value, const int &maxRotationCycles, UpdateFunction &&function)

rotate dimer with trial angle

Rotate the dimer based on a trial rotation in a determined plane and then determine angle of minimum curvature requires 2 calculations per rotation and might not converge fast, but is able to find global minimum more reliable

Return

double Returns the curvature estimate of the rotated dimer

Parameters
  • parameters: The parameters to be optimized

  • value: The value of the optimized point

  • maxRotationCycles: The number of maximum allowed rotations

Template Parameters
  • UpdateFunction: A lambda function with a void return value, and the arguments:

    1. const Eigen::VectorXd& parameters

    2. double& value

    3. Eigen::VectorXd& gradients

template<class UpdateFunction>
double rotationWithGradient(const Eigen::VectorXd &parameters, const double &value, const int &maxRotationCycles, UpdateFunction &&function)

Rotate the dimer directly based on the difference of the forces orthogonal to the dimer axis requires 1 calculation per rotation, only finds closest local minimum.

Return

double Returns the curvature estimate of the rotated dimer

Parameters
  • parameters: The parameters to be optimized

  • value: The value of the optimized point

  • maxRotationCycles: The number of maximum allowed rotations

Template Parameters
  • UpdateFunction: A lambda function with a void return value, and the arguments:

    1. const Eigen::VectorXd& parameters

    2. double& value

    3. Eigen::VectorXd& gradients

template<class UpdateFunction>
double scaleWithProjectionTrustRadius(Eigen::VectorXd &parameters, double &value, UpdateFunction &&function, bool &previousMaxScaling)

Scale the step size based on projection change of the modified force onto the stepvector.

Return

double Returns the curvature estimate of the rotated dimer

Parameters
  • parameters: The parameters to be optimized

  • value: The value of the optimized point

Template Parameters
  • UpdateFunction: A lambda function with a void return value, and the arguments:

    1. const Eigen::VectorXd& parameters

    2. double& value

    3. Eigen::VectorXd& gradients

Parameters
  • previousMaxScaling: The boolean whether maxscaling was used in the last scaling

Private Members

Eigen::VectorXd _dimerAxis
Eigen::VectorXd _gradients
Eigen::VectorXd _parametersR1
Eigen::VectorXd _oldParametersR1
Eigen::VectorXd _gradientsR1
Eigen::VectorXd _fPara
Eigen::VectorXd _fParaOld
Eigen::MatrixXd _dx
Eigen::MatrixXd _dg
Eigen::VectorXd _modForce
Eigen::VectorXd _G
Eigen::VectorXd _orthoG
Eigen::VectorXd _previousOrthoG
Eigen::VectorXd _orthoFN
Eigen::VectorXd _previousOrthoFN
Eigen::VectorXd _steps
unsigned int _numberOfPerformedRotationCycles = 0
unsigned int _sumRotations = 0