IMP  2.0.0
The Integrative Modeling Platform
IMP::atom::MolecularDynamics Class Reference

Simple molecular dynamics optimizer. More...

#include <IMP/atom/MolecularDynamics.h>

+ Inheritance diagram for IMP::atom::MolecularDynamics:

Public Member Functions

 MolecularDynamics (Model *m=nullptr)
 
virtual void assign_velocities (Float temperature)
 Assign velocities representative of the given temperature.
 
virtual Float get_kinetic_energy () const
 
Float get_kinetic_temperature (Float ekinetic) const
 
void set_velocity_cap (Float velocity_cap)
 Set maximum velocity in A/fs. More...
 
- Public Member Functions inherited from IMP::atom::Simulator
 Simulator (Model *m, std::string name="Simulator %1%")
 Create the optimizer. More...
 
double get_current_time () const
 
double get_kt () const
 
ParticleIndexes get_simulation_particle_indexes () const
 
ParticlesTemp get_simulation_particles () const
 
double get_temperature () const
 
 IMP_OPTIMIZER (Simulator)
 
void set_temperature (double d)
 
double simulate (double time_in_fs)
 Simulate until the given time in fs.
 
void set_maximum_time_step (double ts)
 
double get_maximum_time_step () const
 
double get_last_time_step () const
 
void remove_particle (Particle *d)
 
void remove_particles (const Particles &d)
 
void set_particles (const Particles &ps)
 
void set_particles_order (const Particles &objs)
 
unsigned int add_particle (Particle *obj)
 
void add_particles (const Particles &objs)
 
void clear_particles ()
 
unsigned int get_number_of_particles () const
 
bool get_has_particles ()
 
Particleget_particle (unsigned int i) const
 
Particles get_particles () const
 
void reserve_particles (unsigned int sz)
 
- Public Member Functions inherited from IMP::kernel::Optimizer
 Optimizer (Model *m, std::string name="Optimizer %1%")
 
double get_last_score () const
 Return the score found in the last evaluate.
 
Modelget_model () const
 Get the model being optimized.
 
ScoringFunctionget_scoring_function () const
 Return the scoring function that is being used.
 
bool get_stop_on_good_score () const
 
double optimize (unsigned int max_steps)
 Optimize the model for up to max_steps iterations. More...
 
void set_model (Model *m)
 Set the model being optimized. More...
 
virtual void set_scoring_function (ScoringFunctionAdaptor sf)
 
void set_stop_on_good_score (bool tf)
 
virtual void show (std::ostream &out=std::cout) const
 Print info about the optimizer state. More...
 
void remove_optimizer_state (OptimizerState *d)
 
void remove_optimizer_states (const OptimizerStates &d)
 
void set_optimizer_states (const OptimizerStates &ps)
 
void set_optimizer_states_order (const OptimizerStates &objs)
 
unsigned int add_optimizer_state (OptimizerState *obj)
 
void add_optimizer_states (const OptimizerStates &objs)
 
void clear_optimizer_states ()
 
unsigned int get_number_of_optimizer_states () const
 
bool get_has_optimizer_states ()
 
OptimizerStateget_optimizer_state (unsigned int i) const
 
OptimizerStates get_optimizer_states () const
 
void reserve_optimizer_states (unsigned int sz)
 
- Public Member Functions inherited from IMP::base::Object
virtual void clear_caches ()
 
virtual IMP::base::VersionInfo get_version_info () const =0
 Get information about the module and version of the object.
 
void set_check_level (CheckLevel l)
 
void set_log_level (LogLevel l)
 Set the logging level used in this object. More...
 
void set_was_used (bool tf) const
 
void show (std::ostream &out=std::cout) const
 
const std::string & get_name () const
 
void set_name (std::string name)
 

Protected Member Functions

void cap_velocity_component (Float &vel)
 Cap a velocity component to the maximum value.
 
void initialize ()
 
virtual void propagate_coordinates (const ParticleIndexes &ps, double step_size)
 First part of velocity Verlet (update coordinates and half-step velocity)
 
virtual void propagate_velocities (const ParticleIndexes &ps, double step_size)
 Second part of velocity Verlet (update velocity)
 
virtual void setup_degrees_of_freedom (const ParticleIndexes &ps)
 

Protected Attributes

int degrees_of_freedom_
 Number of degrees of freedom in the system.
 
Float velocity_cap_
 Maximum absolute value of a single velocity component.
 
FloatKey vs_ [3]
 Keys of the xyz velocities.
 

Additional Inherited Members

Detailed Description

The particles to be optimized must have optimizable x,y,z attributes and a non-optimizable mass attribute; this optimizer assumes the score to be energy in kcal/mol, the xyz coordinates to be in angstroms, and the mass to be in AMU (g/mol).

Note
RigidBody particles are not handled properly.

Particles without optimized x,y,z and nonoptimized mass are skipped.

See Also
VelocityScalingOptimizerState
LangevinThermostatOptimizerState
BerendsenThermostatOptimizerState
RemoveRigidMotionOptimizerState

Definition at line 34 of file atom/MolecularDynamics.h.

Constructor & Destructor Documentation

IMP::atom::MolecularDynamics::MolecularDynamics ( Model m = nullptr)

Score based on the provided model

Member Function Documentation

virtual Float IMP::atom::MolecularDynamics::get_kinetic_energy ( ) const
virtual
Returns
the current kinetic energy of the system, in kcal/mol

Reimplemented in IMP::isd::MolecularDynamics.

Float IMP::atom::MolecularDynamics::get_kinetic_temperature ( Float  ekinetic) const
Returns
the current kinetic temperature of the system
Parameters
[in]ekinetickinetic energy, e.g. from get_kinetic_energy()
void IMP::atom::MolecularDynamics::set_velocity_cap ( Float  velocity_cap)

At each dynamics time step, the absolute value of each velocity component is capped at this value. This prevents spurious strong forces (occasionally encountered with frustrated conformations) from causing large oscillations in the system. By default, velocities are not capped.

Note
The actual velocities that are capped are the half-step velocities in the velocity Verlet algorithm.

Definition at line 58 of file atom/MolecularDynamics.h.


The documentation for this class was generated from the following file: