OpenMM
 All Classes Namespaces Functions Variables Enumerations Enumerator Friends
ContextImpl Class Reference

This is the internal implementation of a Context. More...

#include <ContextImpl.h>

Public Member Functions

 ContextImpl (Context &owner, System &system, Integrator &integrator, Platform *platform, const std::map< std::string, std::string > &properties)
 Create an ContextImpl for a Context;.
 
 ~ContextImpl ()
 
ContextgetOwner ()
 Get the Context for which this is the implementation.
 
SystemgetSystem ()
 Get System being simulated in this context.
 
IntegratorgetIntegrator ()
 Get Integrator being used to by this context.
 
PlatformgetPlatform ()
 Get the Platform implementation being used for computations.
 
double getTime () const
 Get the current time (in picoseconds).
 
void setTime (double t)
 Set the current time (in picoseconds).
 
void getPositions (std::vector< Vec3 > &positions)
 Get the positions of all particles.
 
void setPositions (const std::vector< Vec3 > &positions)
 Set the positions of all particles.
 
void getVelocities (std::vector< Vec3 > &velocities)
 Get the velocities of all particles.
 
void setVelocities (const std::vector< Vec3 > &velocities)
 Set the velocities of all particles.
 
void getForces (std::vector< Vec3 > &forces)
 Get the current forces on all particles.
 
const std::map< std::string,
double > & 
getParameters () const
 Get the set of all adjustable parameters and their values.
 
double getParameter (std::string name)
 Get the value of an adjustable parameter.
 
void setParameter (std::string name, double value)
 Set the value of an adjustable parameter.
 
void getPeriodicBoxVectors (Vec3 &a, Vec3 &b, Vec3 &c)
 Get the vectors defining the axes of the periodic box (measured in nm).
 
void setPeriodicBoxVectors (const Vec3 &a, const Vec3 &b, const Vec3 &c)
 Set the vectors defining the axes of the periodic box (measured in nm).
 
void applyConstraints (double tol)
 Update the positions of particles so that all distance constraints are satisfied.
 
void applyVelocityConstraints (double tol)
 Update the velocities of particles so the net velocity of each constrained distance is zero.
 
void computeVirtualSites ()
 Recompute the locations of all virtual sites.
 
double calcForcesAndEnergy (bool includeForces, bool includeEnergy, int groups=0xFFFFFFFF)
 Recalculate all of the forces in the system and/or the potential energy of the system (in kJ/mol).
 
int getLastForceGroups () const
 Get the set of force group flags that were passed to the most recent call to calcForcesAndEnergy().
 
double calcKineticEnergy ()
 Calculate the kinetic energy of the system (in kJ/mol).
 
void updateContextState ()
 This should be called at the start of each time step.
 
const std::vector< ForceImpl * > & getForceImpls () const
 Get the list of ForceImpls belonging to this ContextImpl.
 
std::vector< ForceImpl * > & getForceImpls ()
 Get the list of ForceImpls belonging to this ContextImpl.
 
void * getPlatformData ()
 Get the platform-specific data stored in this context.
 
const void * getPlatformData () const
 Get the platform-specific data stored in this context.
 
void setPlatformData (void *data)
 Set the platform-specific data stored in this context.
 
const std::vector< std::vector
< int > > & 
getMolecules () const
 Get a list of the particles in each molecules in the system.
 
void createCheckpoint (std::ostream &stream)
 Create a checkpoint recording the current state of the Context.
 
void loadCheckpoint (std::istream &stream)
 Load a checkpoint that was written by createCheckpoint().
 
void integratorDeleted ()
 This is invoked by the Integrator when it is deleted.
 

Friends

class Context
 

Detailed Description

This is the internal implementation of a Context.

Constructor & Destructor Documentation

ContextImpl ( Context owner,
System system,
Integrator integrator,
Platform platform,
const std::map< std::string, std::string > &  properties 
)

Create an ContextImpl for a Context;.

Member Function Documentation

void applyConstraints ( double  tol)

Update the positions of particles so that all distance constraints are satisfied.

This also recomputes the locations of all virtual sites.

Parameters
tolthe distance tolerance within which constraints must be satisfied.
void applyVelocityConstraints ( double  tol)

Update the velocities of particles so the net velocity of each constrained distance is zero.

Parameters
tolthe velocity tolerance within which constraints must be satisfied.
double calcForcesAndEnergy ( bool  includeForces,
bool  includeEnergy,
int  groups = 0xFFFFFFFF 
)

Recalculate all of the forces in the system and/or the potential energy of the system (in kJ/mol).

After calling this, use getForces() to retrieve the forces that were calculated.

Parameters
includeForcestrue if forces should be calculated
includeEnergytrue if the energy should be calculated
groupsa set of bit flags for which force groups to include. Group i will be included if (groups&(1<<i)) != 0. The default value includes all groups.
Returns
the potential energy of the system, or 0 if includeEnergy is false
double calcKineticEnergy ( )

Calculate the kinetic energy of the system (in kJ/mol).

void computeVirtualSites ( )

Recompute the locations of all virtual sites.

There is rarely a reason to call this, since virtual sites are also updated by applyConstraints(). This is only for the rare situations when you want to enforce virtual sites but not constraints.

void createCheckpoint ( std::ostream &  stream)

Create a checkpoint recording the current state of the Context.

Parameters
streaman output stream the checkpoint data should be written to
const std::vector<ForceImpl*>& getForceImpls ( ) const

Get the list of ForceImpls belonging to this ContextImpl.

std::vector<ForceImpl*>& getForceImpls ( )

Get the list of ForceImpls belonging to this ContextImpl.

void getForces ( std::vector< Vec3 > &  forces)

Get the current forces on all particles.

Parameters
forceson exit, this contains the forces
Integrator& getIntegrator ( )
inline

Get Integrator being used to by this context.

int getLastForceGroups ( ) const

Get the set of force group flags that were passed to the most recent call to calcForcesAndEnergy().

const std::vector<std::vector<int> >& getMolecules ( ) const

Get a list of the particles in each molecules in the system.

Two particles are in the same molecule if they are connected by constraints or bonds.

Context& getOwner ( )
inline

Get the Context for which this is the implementation.

double getParameter ( std::string  name)

Get the value of an adjustable parameter.

If there is no parameter with the specified name, this throws an exception.

Parameters
namethe name of the parameter to get
const std::map<std::string, double>& getParameters ( ) const

Get the set of all adjustable parameters and their values.

void getPeriodicBoxVectors ( Vec3 a,
Vec3 b,
Vec3 c 
)

Get the vectors defining the axes of the periodic box (measured in nm).

They will affect any Force that uses periodic boundary conditions.

Currently, only rectangular boxes are supported. This means that a, b, and c must be aligned with the x, y, and z axes respectively. Future releases may support arbitrary triclinic boxes.

Parameters
athe vector defining the first edge of the periodic box
bthe vector defining the second edge of the periodic box
cthe vector defining the third edge of the periodic box
Platform& getPlatform ( )
inline

Get the Platform implementation being used for computations.

void* getPlatformData ( )

Get the platform-specific data stored in this context.

const void* getPlatformData ( ) const

Get the platform-specific data stored in this context.

void getPositions ( std::vector< Vec3 > &  positions)

Get the positions of all particles.

Parameters
positionson exit, this contains the particle positions
System& getSystem ( )
inline

Get System being simulated in this context.

double getTime ( ) const

Get the current time (in picoseconds).

void getVelocities ( std::vector< Vec3 > &  velocities)

Get the velocities of all particles.

Parameters
velocitieson exit, this contains the particle velocities
void integratorDeleted ( )
inline

This is invoked by the Integrator when it is deleted.

This is needed to ensure the cleanup process is done correctly, since we don't know whether the Integrator or Context will be deleted first.

void loadCheckpoint ( std::istream &  stream)

Load a checkpoint that was written by createCheckpoint().

Parameters
streaman input stream the checkpoint data should be read from
void setParameter ( std::string  name,
double  value 
)

Set the value of an adjustable parameter.

If there is no parameter with the specified name, this throws an exception.

Parameters
namethe name of the parameter to set
valuethe value of the parameter
void setPeriodicBoxVectors ( const Vec3 a,
const Vec3 b,
const Vec3 c 
)

Set the vectors defining the axes of the periodic box (measured in nm).

They will affect any Force that uses periodic boundary conditions.

Currently, only rectangular boxes are supported. This means that a, b, and c must be aligned with the x, y, and z axes respectively. Future releases may support arbitrary triclinic boxes.

Parameters
athe vector defining the first edge of the periodic box
bthe vector defining the second edge of the periodic box
cthe vector defining the third edge of the periodic box
void setPlatformData ( void *  data)

Set the platform-specific data stored in this context.

void setPositions ( const std::vector< Vec3 > &  positions)

Set the positions of all particles.

Parameters
positionsa vector containg the particle positions
void setTime ( double  t)

Set the current time (in picoseconds).

void setVelocities ( const std::vector< Vec3 > &  velocities)

Set the velocities of all particles.

Parameters
velocitiesa vector containg the particle velocities
void updateContextState ( )

This should be called at the start of each time step.

It calls updateContextState() on each ForceImpl in the system, allowing them to modify the values of state variables.

Friends And Related Function Documentation

friend class Context
friend

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