class atmd::AtomicDimer¶
Overview¶
Atomic dimer. More…
#include <AtomicDimer.h>
class AtomicDimer {
public:
// structs
struct VarianceGateConfig;
struct early_stopping_t;
// fields
gpr::AtomsConfiguration atoms_config;
// construction
AtomicDimer();
virtual ~AtomicDimer();
// methods
void initialize(const gpr::InputParameters& parameters, const gpr::Observation& initial_data, const gpr::Observation& middle_point, const gpr::Coord& initial_orientation, const gpr::AtomsConfiguration& init_atom_config);
void setAtomsConfiguration(const gpr::AtomsConfiguration& _atoms_config);
void execute(pot::PotentialWrapper& general_potential);
Eigen::VectorXd getFinalForceAtMidPoint() const;
Eigen::VectorXd getFinalCoordOfMidPoint() const;
gpr::Index_t getTotalForceCalls() const;
gpr::Index_t getTotalGPRForceCalls() const;
gpr::Index_t getIterations() const;
double getFinalEnergy() const;
double getFinalCurvature() const;
const gpr::Coord& getFinalOrientation() const;
gpr::Coord& getMidPointCoordForAllObservationPointsRef();
gpr::FieldMatrixd& getEnergyForAllObservationPointsRef();
gpr::Coord& getGradientForAllObservationPointsRef();
gpr::GaussianProcessRegression* getGPRModel();
void callGeneralPotentialFromEON(const gpr::vector3_reg(&) cell_dimensions[3], gpr::Observation& mid_point, pot::PotentialWrapper& potential);
bool isFinalConvergenceReached(const gpr::Observation& middle_point, const double stop_criteria);
bool isRotationalConvergenceReached(const gpr::Coord& orient, gpr::Observation& middle_point);
bool rotateDimer(const double stop_criteria, gpr::Coord& orient, gpr::LBFGSInfo& rot_info, gpr::Observation& middle_point, gpr::Observation& approx_obs);
void translateDimer(const gpr::Observation& middle_point, const gpr::Coord& orient, const gpr::Observation& approx_obs, const double curvature, gpr::LBFGSInfo& trans_info, gpr::Coord& R_new);
bool updateActiveFrozenAtoms(const gpr::Coord& R_new, gpr::LBFGSInfo& rot_info, gpr::LBFGSInfo& trans_info);
void setVarianceGateConfig(const VarianceGateConfig& cfg);
const VarianceGateConfig& getVarianceGateConfig() const;
void setVarianceAcceptTau(double tau);
double getVarianceAcceptTau() const;
void setVarianceAcceptTauRelative(double tau_relative);
double getVarianceAcceptTauRelative() const;
gpr::Index_t getVarianceSkippedAccurateCalls() const;
};
Detailed Documentation¶
Atomic dimer.
Fields¶
gpr::AtomsConfiguration atoms_config
Configuration of atoms
Construction¶
AtomicDimer()
Default constructor.
Methods¶
void initialize(const gpr::InputParameters& parameters, const gpr::Observation& initial_data, const gpr::Observation& middle_point, const gpr::Coord& initial_orientation, const gpr::AtomsConfiguration& init_atom_config)
Initilalize with initial data.
Parameters:
parameters |
Structure of parameters from the input file |
initial_data |
Initial observation (energy, gradients, positions) |
middle_point |
Initial observation at the middle point (energy, gradients, positions) |
initial_orientation |
Initial orientation of the dimer |
init_atom_config |
Structure with initial configuration of atoms |
void setAtomsConfiguration(const gpr::AtomsConfiguration& _atoms_config)
Set configuration of the atoms.
Parameters:
_atoms_config |
A reference configuration |
void execute(pot::PotentialWrapper& general_potential)
This function uses the atomic GP-dimer method to converge to a saddle point, starting from somewhere inside the convergence area.
The relaxation of the dimer on the approximated energy surface is done according to a dimer method, where a rotation step rotates the dimer (a pair of images) towards its minimum energy orientation to find the lowest curvature mode of the potential energy and translation step moves the dimer towards the saddle point by inverting the force component in the direction of the dimer. After each relaxation phase, the energy and gradient are acquired at the middle point of the dimer, and the GP hyperparameters are reoptimized.
Parameters:
general_potential |
Reference to a class of general potential or to the object of Matter class (in a case of coupling with EON) |
Eigen::VectorXd getFinalForceAtMidPoint() const
Return final force at the middle point.
Eigen::VectorXd getFinalCoordOfMidPoint() const
Return final coordinates of the middle point.
gpr::Index_t getTotalForceCalls() const
Return the total number of force calls.
gpr::Index_t getTotalGPRForceCalls() const
Return the total number of force calls.
gpr::Index_t getIterations() const
Return the total number of cycles.
This is essentially the number of relaxation phases
double getFinalEnergy() const
Return final energy.
double getFinalCurvature() const
Return final curvature.
const gpr::Coord& getFinalOrientation() const
Return final orientation of the dimer.
gpr::Coord& getMidPointCoordForAllObservationPointsRef()
Return gpr::Coordinates of all observed points.
gpr::FieldMatrixd& getEnergyForAllObservationPointsRef()
Return energies for all observed points.
Return energies for all observation points.
gpr::Coord& getGradientForAllObservationPointsRef()
Return gradients for all observed points.
Return gradients for all observation points.
gpr::GaussianProcessRegression* getGPRModel()
Get GPR model.
void callGeneralPotentialFromEON(const gpr::vector3_reg(&) cell_dimensions[3], gpr::Observation& mid_point, pot::PotentialWrapper& potential)
Call for a potential defined in a Matter object.
Parameters:
cell_dimensions |
Dimensions of the cell in EON 9 element format |
mid_point |
The point on which to obtain the potential |
T |
The potential object pointer |
bool isFinalConvergenceReached(const gpr::Observation& middle_point, const double stop_criteria)
Final convergence test for force.
@detail Function find maximum absolute component of the force acting on the middle point of the dimer and returns true if the value is less than T_dimer.
Parameters:
middle_point |
The current mid-point of the configuration |
stop_criteria |
The final convergence threshold |
Returns:
true if final convergence is reached
bool isRotationalConvergenceReached(const gpr::Coord& orient, gpr::Observation& middle_point)
Rotational convergence test.
Parameters:
orient |
The orientation coordinates |
stop_criteria |
The final convergence threshold |
Returns:
true if rotational convergence is reached
bool rotateDimer(const double stop_criteria, gpr::Coord& orient, gpr::LBFGSInfo& rot_info, gpr::Observation& middle_point, gpr::Observation& approx_obs)
Rotate the dimer.
void translateDimer(const gpr::Observation& middle_point, const gpr::Coord& orient, const gpr::Observation& approx_obs, const double curvature, gpr::LBFGSInfo& trans_info, gpr::Coord& R_new)
Translate the dimer.
bool updateActiveFrozenAtoms(const gpr::Coord& R_new, gpr::LBFGSInfo& rot_info, gpr::LBFGSInfo& trans_info)
Update active frozen atoms.
Returns:
True if new atoms are activated
void setVarianceAcceptTau(double tau)
Convenience setter for the absolute accept_tau alone; other fields of the default-constructed VarianceGateConfig are preserved.
void setVarianceAcceptTauRelative(double tau_relative)
Convenience setter for the relative accept_tau alone. The effective absolute threshold is recomputed on every gate call from the current magnSigma2, so this setter need not be re-invoked when SCG re-fits the hyperparameters between outer iterations.