|
double | get_w () |
|
int | get_n_hist_max () |
|
double | get_initial_step_size () |
|
| periodic_optimizer (int &nat) |
| Construct a new periodic optimizer::periodic optimizer object for free or fixed cell optimization with default parameters. More...
|
|
| periodic_optimizer (int &nat, double initial_step_size, int nhist_max, double alpha0, double eps_subsp) |
| Construct a new periodic optimizer::periodic optimizer object for free or fixed cell optimization with custom parameters. More...
|
|
| periodic_optimizer (int &nat, Eigen::Vector3d &lat_a, Eigen::Vector3d &lat_b, Eigen::Vector3d &lat_c) |
| Construct a new periodic optimizer::periodic optimizer object for variable cell shape optimization with default parameters. More...
|
|
| periodic_optimizer (int &nat, Eigen::Vector3d &lat_a, Eigen::Vector3d &lat_b, Eigen::Vector3d &lat_c, double initial_step_size, int nhist_max, double lattice_weight, double alpha0, double eps_subsp) |
| Construct a new periodic optimizer::periodic optimizer object for variable cell shape optimization with custom parameters. More...
|
|
void | step (Eigen::MatrixXd &r, double &energy, Eigen::MatrixXd &f) |
| Calculates new atomic coordinates that are closer to the local minimum. Fixed cell optimization. This function should be used the following way: More...
|
|
void | step (Eigen::MatrixXd &r, double &energy, Eigen::MatrixXd &f, Eigen::Vector3d &lat_a, Eigen::Vector3d &lat_b, Eigen::Vector3d &lat_c, Eigen::Matrix3d &stress) |
| Calculates new atomic coordinates that are closer to the local minimum. Variable cell shape optimization. This function should be used the following way: More...
|
|
void | check_forces (Eigen::MatrixXd forces) |
|
double | lower_bound () |
| Estimates a lower bound of the energy of the local minimum. More...
|
|
Definition at line 20 of file periodic_optimizer.hpp.
◆ periodic_optimizer() [1/4]
vcsqnm::PES_optimizer::periodic_optimizer::periodic_optimizer |
( |
int & |
nat | ) |
|
|
inline |
Construct a new periodic optimizer::periodic optimizer object for free or fixed cell optimization with default parameters.
- Parameters
-
Definition at line 55 of file periodic_optimizer.hpp.
◆ periodic_optimizer() [2/4]
vcsqnm::PES_optimizer::periodic_optimizer::periodic_optimizer |
( |
int & |
nat, |
|
|
double |
initial_step_size, |
|
|
int |
nhist_max, |
|
|
double |
alpha0, |
|
|
double |
eps_subsp |
|
) |
| |
|
inline |
Construct a new periodic optimizer::periodic optimizer object for free or fixed cell optimization with custom parameters.
- Parameters
-
nat | number of atoms |
initial_step_size | 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. If a system only contains weaker bonds a value up to 5.0 may speed up the convergence. |
nhist_max | Maximal number of steps that will be stored in the history list. Use a value between 3 and 20. Must be <= than 3*nat. |
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 73 of file periodic_optimizer.hpp.
◆ periodic_optimizer() [3/4]
vcsqnm::PES_optimizer::periodic_optimizer::periodic_optimizer |
( |
int & |
nat, |
|
|
Eigen::Vector3d & |
lat_a, |
|
|
Eigen::Vector3d & |
lat_b, |
|
|
Eigen::Vector3d & |
lat_c |
|
) |
| |
|
inline |
Construct a new periodic optimizer::periodic optimizer object for variable cell shape optimization with default parameters.
- Parameters
-
nat | number of atoms |
lat_a | first lattice vector |
lat_b | second lattice vector |
lat_c | third lattice vector |
Definition at line 91 of file periodic_optimizer.hpp.
◆ periodic_optimizer() [4/4]
vcsqnm::PES_optimizer::periodic_optimizer::periodic_optimizer |
( |
int & |
nat, |
|
|
Eigen::Vector3d & |
lat_a, |
|
|
Eigen::Vector3d & |
lat_b, |
|
|
Eigen::Vector3d & |
lat_c, |
|
|
double |
initial_step_size, |
|
|
int |
nhist_max, |
|
|
double |
lattice_weight, |
|
|
double |
alpha0, |
|
|
double |
eps_subsp |
|
) |
| |
|
inline |
Construct a new periodic optimizer::periodic optimizer object for variable cell shape optimization with custom parameters.
- Parameters
-
nat | number of atoms |
lat_a | first lattice vector |
lat_b | second lattice vector |
lat_c | third lattice vector |
initial_step_size | 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. If a system only contains weaker bonds a value up to 5.0 may speed up the convergence. |
nhist_max | Maximal number of steps that will be stored in the history list. Use a value between 3 and 20. Must be <= than 3*nat + 9. |
lattice_weight | weight / size of the supercell that is used to transform lattice derivatives. Use a value between 1 and 2. Default is 2. |
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 111 of file periodic_optimizer.hpp.
◆ get_w()
double vcsqnm::PES_optimizer::periodic_optimizer::get_w |
( |
| ) |
|
|
inline |
◆ get_n_hist_max()
int vcsqnm::PES_optimizer::periodic_optimizer::get_n_hist_max |
( |
| ) |
|
|
inline |
◆ get_initial_step_size()
double vcsqnm::PES_optimizer::periodic_optimizer::get_initial_step_size |
( |
| ) |
|
|
inline |
◆ step() [1/2]
void vcsqnm::PES_optimizer::periodic_optimizer::step |
( |
Eigen::MatrixXd & |
r, |
|
|
double & |
energy, |
|
|
Eigen::MatrixXd & |
f |
|
) |
| |
|
inline |
Calculates new atomic coordinates that are closer to the local minimum. Fixed cell optimization. This function should be used the following way:
- calculate energies and forces at positions r.
- call the step function to update positions r.
- repeat.
- Parameters
-
r | Input: atomic coordinates, dimension(3, nat). Output: improved coordinates that are calculated based on forces from this and previous iterations. |
energy | Potential energy of the system in it's current state |
f | Forces of the system in it's current state. dimension(3, nat) |
Definition at line 130 of file periodic_optimizer.hpp.
◆ step() [2/2]
void vcsqnm::PES_optimizer::periodic_optimizer::step |
( |
Eigen::MatrixXd & |
r, |
|
|
double & |
energy, |
|
|
Eigen::MatrixXd & |
f, |
|
|
Eigen::Vector3d & |
lat_a, |
|
|
Eigen::Vector3d & |
lat_b, |
|
|
Eigen::Vector3d & |
lat_c, |
|
|
Eigen::Matrix3d & |
stress |
|
) |
| |
|
inline |
Calculates new atomic coordinates that are closer to the local minimum. Variable cell shape optimization. This function should be used the following way:
- calculate energies, forces and stress tensor at positions r and lattice vectors a, b, c.
- call the step function to update positions r and lattice vectors.
- repeat.
- Parameters
-
r | Input: atomic coordinates, dimension(3, nat). Output: improved coordinates that are calculated based on forces from this and previous iterations. |
energy | Potential energy of the system in it's current state |
f | Forces of the system in it's current state. dimension(3, nat) |
lat_a | first lattice vector |
lat_b | second lattice vector |
lat_c | third lattice vector |
stress | stress tensor of the system in it' current state. |
Definition at line 157 of file periodic_optimizer.hpp.
◆ check_forces()
void vcsqnm::PES_optimizer::periodic_optimizer::check_forces |
( |
Eigen::MatrixXd |
forces | ) |
|
|
inline |
◆ lower_bound()
double vcsqnm::PES_optimizer::periodic_optimizer::lower_bound |
( |
| ) |
|
|
inline |
Estimates a lower bound of the energy of the local minimum.
- Returns
- double Lower bound estimate
Definition at line 213 of file periodic_optimizer.hpp.
◆ combine_matrices()
Eigen::VectorXd vcsqnm::PES_optimizer::periodic_optimizer::combine_matrices |
( |
Eigen::MatrixXd |
a, |
|
|
Eigen::MatrixXd |
b |
|
) |
| |
|
inlineprivate |
◆ split_matrices()
void vcsqnm::PES_optimizer::periodic_optimizer::split_matrices |
( |
Eigen::MatrixXd & |
a, |
|
|
Eigen::MatrixXd & |
b, |
|
|
Eigen::VectorXd & |
c |
|
) |
| |
|
inlineprivate |
◆ calc_lattice_derivatices()
Eigen::Matrix3d vcsqnm::PES_optimizer::periodic_optimizer::calc_lattice_derivatices |
( |
Eigen::Matrix3d & |
stress, |
|
|
Eigen::Matrix3d & |
alat |
|
) |
| |
|
inlineprivate |
◆ setupInitialLattice()
void vcsqnm::PES_optimizer::periodic_optimizer::setupInitialLattice |
( |
int & |
nat, |
|
|
Eigen::Vector3d & |
lat_a, |
|
|
Eigen::Vector3d & |
lat_b, |
|
|
Eigen::Vector3d & |
lat_c |
|
) |
| |
|
inlineprivate |
◆ initial_latttice
Eigen::Matrix3d vcsqnm::PES_optimizer::periodic_optimizer::initial_latttice |
|
private |
◆ initial_latttice_inv
Eigen::Matrix3d vcsqnm::PES_optimizer::periodic_optimizer::initial_latttice_inv |
|
private |
◆ lattice_transformer
Eigen::Matrix3d vcsqnm::PES_optimizer::periodic_optimizer::lattice_transformer |
|
private |
◆ lattice_transformer_inv
Eigen::Matrix3d vcsqnm::PES_optimizer::periodic_optimizer::lattice_transformer_inv |
|
private |
◆ opt
std::unique_ptr<sqnm_space::SQNM> vcsqnm::PES_optimizer::periodic_optimizer::opt |
|
private |
◆ nat
int vcsqnm::PES_optimizer::periodic_optimizer::nat |
|
private |
◆ ndim
int vcsqnm::PES_optimizer::periodic_optimizer::ndim |
|
private |
◆ opt_lattice
bool vcsqnm::PES_optimizer::periodic_optimizer::opt_lattice |
|
private |
◆ initial_step_size
double vcsqnm::PES_optimizer::periodic_optimizer::initial_step_size = 1 |
|
private |
◆ n_hist_max
int vcsqnm::PES_optimizer::periodic_optimizer::n_hist_max = 10 |
|
private |
double vcsqnm::PES_optimizer::periodic_optimizer::w = 2.0 |
|
private |
◆ f_std_deviation
double vcsqnm::PES_optimizer::periodic_optimizer::f_std_deviation = 0.0 |
|
private |
The documentation for this class was generated from the following file: