|
| SQNM (int ndim_, int nhistx_, double alpha_) |
| Construct a new SQNM::SQNM object using default parameters. More...
|
|
| SQNM (int ndim_, int nhistx_, double alpha_, double alpha0_, double eps_subsp_) |
| Construct a new SQNM::SQNM object using custom parameters. More...
|
|
Eigen::VectorXd | step (Eigen::VectorXd &x, double &f_of_x, Eigen::VectorXd &df_dx) |
| Calculates new coordinates that are closer to local minimum that the current coordinates. This function should be used the following way: More...
|
|
double | lower_bound () |
| Estimates a lower bound of the energy of the local minimum. More...
|
|
Definition at line 20 of file sqnm.hpp.
◆ SQNM() [1/2]
vcsqnm::sqnm_space::SQNM::SQNM |
( |
int |
ndim_, |
|
|
int |
nhistx_, |
|
|
double |
alpha_ |
|
) |
| |
|
inline |
Construct a new SQNM::SQNM object using default parameters.
- Parameters
-
ndim_ | number of dimensions of target function |
nhistx_ | Maximal number of steps that will be stored in the history list. Use a value between 3 and 20. Must be <= than ndim_. |
alpha_ | initial step size. default is 1.0. For systems with hard bonds (e.g. C-C) use a value between and 1.0 and 2.5 Should be approximately the inverse of the largest eigenvalue of the Hessian matrix. If alpha is negative, the inial step size is estimated using the mechanism from section 6.4 of the vc-sqnm paper: https://arxiv.org/abs/2206.07339 beta will then be equal to minus alpha. Good choices for beta are 0.1 in hartee / bohr^2 and 0.001 in eV / A^2 |
Definition at line 56 of file sqnm.hpp.
◆ SQNM() [2/2]
vcsqnm::sqnm_space::SQNM::SQNM |
( |
int |
ndim_, |
|
|
int |
nhistx_, |
|
|
double |
alpha_, |
|
|
double |
alpha0_, |
|
|
double |
eps_subsp_ |
|
) |
| |
|
inline |
Construct a new SQNM::SQNM object using custom parameters.
- Parameters
-
ndim_ | number of dimensions of target function |
nhistx_ | Maximal number of steps that will be stored in the history list. Use a value between 3 and 20. Must be <= than ndim_. |
alpha_ | initial step size. default is 1.0. For systems with hard bonds (e.g. C-C) use a value between and 1.0 and 2.5. Should be approximately the inverse of the largest eigenvalue of the Hessian matrix. If alpha is negative, the inial step size is estimated using the mechanism from section 6.4 of the vc-sqnm paper: https://arxiv.org/abs/2206.07339 beta will then be equal to minus alpha. Good choices for beta are 0.1 in hartee / bohr^2 and 0.001 in eV / A^2 |
alpha0_ | * |
alpha0 | Lower limit on the step size. 1.e-2 is the default. |
eps_subsp_ | Lower limit on linear dependencies of basis vectors in history list. Default 1.e-4. |
Definition at line 86 of file sqnm.hpp.
◆ step()
Eigen::VectorXd vcsqnm::sqnm_space::SQNM::step |
( |
Eigen::VectorXd & |
x, |
|
|
double & |
f_of_x, |
|
|
Eigen::VectorXd & |
df_dx |
|
) |
| |
|
inline |
Calculates new coordinates that are closer to local minimum that the current coordinates. This function should be used the following way:
- calculate f(x) and the derivative.
- call the step function.
- add return value of step function to x.
- repeat.
- Parameters
-
x | Current position vector |
f_of_x | value of the target function evaluated at position x. |
df_dx | derivative of the target function evaluated at x. |
- Returns
- VectorXd displacent that can be added to x in order to get new improved coordinates.
Definition at line 115 of file sqnm.hpp.
◆ lower_bound()
double vcsqnm::sqnm_space::SQNM::lower_bound |
( |
| ) |
|
|
inline |
Estimates a lower bound of the energy of the local minimum.
- Returns
- double Lower bound estimate
Definition at line 245 of file sqnm.hpp.
◆ calc_gainratio()
double vcsqnm::sqnm_space::SQNM::calc_gainratio |
( |
double & |
f | ) |
|
|
inlineprivate |
◆ adjust_stepsize()
void vcsqnm::sqnm_space::SQNM::adjust_stepsize |
( |
double & |
gainratio | ) |
|
|
inlineprivate |
◆ calc_ovrlp()
Eigen::MatrixXd vcsqnm::sqnm_space::SQNM::calc_ovrlp |
( |
| ) |
|
|
inlineprivate |
◆ ndim
int vcsqnm::sqnm_space::SQNM::ndim |
|
private |
◆ nhistx
int vcsqnm::sqnm_space::SQNM::nhistx |
|
private |
◆ eps_subsp
double vcsqnm::sqnm_space::SQNM::eps_subsp = 1.e-3 |
|
private |
◆ alpha0
double vcsqnm::sqnm_space::SQNM::alpha0 = 1.e-2 |
|
private |
◆ xlist
◆ flist
◆ alpha
double vcsqnm::sqnm_space::SQNM::alpha |
|
private |
◆ dir_of_descent
Eigen::VectorXd vcsqnm::sqnm_space::SQNM::dir_of_descent |
|
private |
◆ prev_f
double vcsqnm::sqnm_space::SQNM::prev_f |
|
private |
◆ prev_df_dx
Eigen::VectorXd vcsqnm::sqnm_space::SQNM::prev_df_dx |
|
private |
◆ expected_positions
Eigen::VectorXd vcsqnm::sqnm_space::SQNM::expected_positions |
|
private |
◆ h_subsp
Eigen::MatrixXd vcsqnm::sqnm_space::SQNM::h_subsp |
|
private |
◆ h_evec_subsp
Eigen::MatrixXd vcsqnm::sqnm_space::SQNM::h_evec_subsp |
|
private |
◆ h_evec
Eigen::MatrixXd vcsqnm::sqnm_space::SQNM::h_evec |
|
private |
◆ h_eval
Eigen::VectorXd vcsqnm::sqnm_space::SQNM::h_eval |
|
private |
◆ res
Eigen::VectorXd vcsqnm::sqnm_space::SQNM::res |
|
private |
◆ esolve
Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> vcsqnm::sqnm_space::SQNM::esolve |
|
private |
◆ res_temp
Eigen::VectorXd vcsqnm::sqnm_space::SQNM::res_temp |
|
private |
◆ nhist
int vcsqnm::sqnm_space::SQNM::nhist = 0 |
|
private |
◆ estimate_step_size
bool vcsqnm::sqnm_space::SQNM::estimate_step_size = false |
|
private |
The documentation for this class was generated from the following file: