1 #ifndef OPENMM_RPMDINTEGRATOR_H_
2 #define OPENMM_RPMDINTEGRATOR_H_
35 #include "openmm/Integrator.h"
36 #include "openmm/Kernel.h"
37 #include "openmm/State.h"
38 #include "openmm/Vec3.h"
39 #include "openmm/internal/windowsExportRpmd.h"
79 RPMDIntegrator(
int numCopies,
double temperature,
double frictionCoeff,
double stepSize);
92 RPMDIntegrator(
int numCopies,
double temperature,
double frictionCoeff,
double stepSize,
const std::map<int, int>& contractions);
137 return applyThermostat;
143 applyThermostat = apply;
149 return randomNumberSeed;
164 randomNumberSeed = seed;
181 void setPositions(
int copy,
const std::vector<Vec3>& positions);
188 void setVelocities(
int copy,
const std::vector<Vec3>& velocities);
201 State getState(
int copy,
int types,
bool enforcePeriodicBox=
false,
int groups=0xFFFFFFFF);
206 double getTotalEnergy();
212 void step(
int steps);
219 void initialize(ContextImpl& context);
232 std::vector<std::string> getKernelNames();
236 double computeKineticEnergy();
238 double temperature, friction;
239 int numCopies, randomNumberSeed;
240 bool applyThermostat;
241 std::map<int, int> contractions;
242 bool forcesAreValid, hasSetPosition, hasSetVelocity, isFirstStep;
const std::map< int, int > & getContractions() const
Get the ring polymer contractions to use for evaluating different force groups.
Definition: RPMDIntegrator.h:172
void setFriction(double coeff)
Set the friction coefficient which determines how strongly the system is coupled to the heat bath (in...
Definition: RPMDIntegrator.h:130
void setTemperature(double temp)
Set the temperature of the heat bath (in Kelvin).
Definition: RPMDIntegrator.h:112
bool getApplyThermostat() const
Get whether a thermostat is applied to the system.
Definition: RPMDIntegrator.h:136
A State object records a snapshot of the current state of a simulation at a point in time...
Definition: State.h:54
int getRandomNumberSeed() const
Get the random number seed.
Definition: RPMDIntegrator.h:148
double getTemperature() const
Get the temperature of the heat bath (in Kelvin).
Definition: RPMDIntegrator.h:104
int getNumCopies() const
Get the number of copies of the system being simulated.
Definition: RPMDIntegrator.h:96
An Integrator defines a method for simulating a System by integrating the equations of motion...
Definition: Integrator.h:54
This is an Integrator which simulates a System using ring polymer molecular dynamics (RPMD)...
Definition: RPMDIntegrator.h:69
DataType
This is an enumeration of the types of data which may be stored in a State.
Definition: State.h:61
A Kernel encapsulates a particular implementation of a calculation that can be performed on the data ...
Definition: Kernel.h:58
double getFriction() const
Get the friction coefficient which determines how strongly the system is coupled to the heat bath (in...
Definition: RPMDIntegrator.h:121
Definition: AndersenThermostat.h:40
void setRandomNumberSeed(int seed)
Set the random number seed.
Definition: RPMDIntegrator.h:163
void setApplyThermostat(bool apply)
Set whether a thermostat is applied to the system.
Definition: RPMDIntegrator.h:142