OpenMM
 All Classes Namespaces Functions Variables Pages
System Class Reference

This class represents a molecular system. More...

Inherits _object.

Public Member Functions

def __del__
 del(OpenMM::System self) More...
 
def getNumParticles
 getNumParticles(System self) -> int More...
 
def addParticle
 addParticle(System self, double mass) -> int More...
 
def getParticleMass
 getParticleMass(System self, int index) -> double More...
 
def setParticleMass
 setParticleMass(System self, int index, double mass) More...
 
def setVirtualSite
 setVirtualSite(System self, int index, VirtualSite virtualSite) More...
 
def isVirtualSite
 isVirtualSite(System self, int index) -> bool More...
 
def getVirtualSite
 getVirtualSite(System self, int index) -> VirtualSite More...
 
def getNumConstraints
 getNumConstraints(System self) -> int More...
 
def addConstraint
 addConstraint(System self, int particle1, int particle2, double distance) -> int More...
 
def getConstraintParameters
 getConstraintParameters(System self, int index) More...
 
def setConstraintParameters
 setConstraintParameters(System self, int index, int particle1, int particle2, double distance) More...
 
def addForce
 addForce(System self, Force force) -> int More...
 
def getNumForces
 getNumForces(System self) -> int More...
 
def getForce
 getForce(System self, int index) -> Force getForce(System self, int index) -> Force More...
 
def getDefaultPeriodicBoxVectors
 getDefaultPeriodicBoxVectors(System self) More...
 
def setDefaultPeriodicBoxVectors
 setDefaultPeriodicBoxVectors(System self, Vec3 const & a, Vec3 const & b, Vec3 const & c) More...
 
def __getstate__
 
def __setstate__
 
def getForces
 Get the list of Forces in this System. More...
 
def __init__
 init(OpenMM::System self) -> System init(OpenMM::System self, System other) -> System More...
 

Public Attributes

 this
 

Detailed Description

This class represents a molecular system.

The definition of a System involves four elements:

The particles and constraints are defined directly by the System object, while forces are defined by objects that extend the Force class. After creating a System, call addParticle() once for each particle, addConstraint() for each constraint, and addForce() for each Force.

In addition, particles may be designated as "virtual sites". These are particles whose positions are computed automatically based on the positions of other particles. To define a virtual site, call setVirtualSite(), passing in a VirtualSite object that defines the rules for computing its position.

Constructor & Destructor Documentation

def __del__ (   self)

del(OpenMM::System self)

References simtk.openmm.openmm.stripUnits().

def __init__ (   self,
  args 
)

init(OpenMM::System self) -> System init(OpenMM::System self, System other) -> System

Create a new System.

References simtk.openmm.openmm.stripUnits(), ios.this, ostream.this, istream.this, iostream.this, pairii.this, vectord.this, vectorddd.this, vectori.this, vectorii.this, vectorpairii.this, vectorstring.this, mapstringstring.this, mapstringdouble.this, mapii.this, seti.this, AmoebaAngleForce.this, AmoebaBondForce.this, AmoebaGeneralizedKirkwoodForce.this, AmoebaInPlaneAngleForce.this, AmoebaMultipoleForce.this, AmoebaOutOfPlaneBendForce.this, AmoebaPiTorsionForce.this, AmoebaStretchBendForce.this, AmoebaTorsionTorsionForce.this, AmoebaVdwForce.this, AmoebaWcaDispersionForce.this, AndersenThermostat.this, BrownianIntegrator.this, CMAPTorsionForce.this, CMMotionRemover.this, Context.this, Continuous1DFunction.this, Continuous2DFunction.this, Continuous3DFunction.this, CustomAngleForce.this, CustomBondForce.this, CustomCompoundBondForce.this, CustomExternalForce.this, CustomGBForce.this, CustomHbondForce.this, CustomIntegrator.this, CustomNonbondedForce.this, CustomTorsionForce.this, Discrete1DFunction.this, Discrete2DFunction.this, Discrete3DFunction.this, DrudeForce.this, DrudeLangevinIntegrator.this, DrudeSCFIntegrator.this, GBSAOBCForce.this, GBVIForce.this, HarmonicAngleForce.this, HarmonicBondForce.this, LangevinIntegrator.this, LocalCoordinatesSite.this, MonteCarloAnisotropicBarostat.this, MonteCarloBarostat.this, NonbondedForce.this, OpenMMException.this, OutOfPlaneSite.this, PeriodicTorsionForce.this, RBTorsionForce.this, RPMDIntegrator.this, and System.this.

Member Function Documentation

def __getstate__ (   self)
def __setstate__ (   self,
  serializationString 
)
def addConstraint (   self,
  args 
)

addConstraint(System self, int particle1, int particle2, double distance) -> int

Add a constraint to the System. Particles whose mass is 0 cannot participate in constraints.

Parameters
particle1the index of the first particle involved in the constraint
particle2the index of the second particle involved in the constraint
distancethe required distance between the two particles, measured in nm

References simtk.openmm.openmm.stripUnits().

def addForce (   self,
  args 
)

addForce(System self, Force force) -> int

Add a Force to the System. The Force should have been created on the heap with the "new" operator. The System takes over ownership of it, and deletes the Force when the System itself is deleted.

Parameters
forcea pointer to the Force object to be added
def addParticle (   self,
  args 
)

addParticle(System self, double mass) -> int

Add a particle to the System. If the mass is 0, Integrators will ignore the particle and not modify its position or velocity. This is most often used for virtual sites, but can also be used as a way to prevent a particle from moving.

Parameters
massthe mass of the particle (in atomic mass units)

References simtk.openmm.openmm.stripUnits().

def getConstraintParameters (   self,
  args 
)

getConstraintParameters(System self, int index)

Get the parameters defining a distance constraint.

Parameters
indexthe index of the constraint for which to get parameters
particle1the index of the first particle involved in the constraint
particle2the index of the second particle involved in the constraint
distancethe required distance between the two particles, measured in nm

References simtk.openmm.openmm.stripUnits().

def getDefaultPeriodicBoxVectors (   self)

getDefaultPeriodicBoxVectors(System self)

Get the default values of the vectors defining the axes of the periodic box (measured in nm). Any newly created Context will have its box vectors set to these. They will affect any Force added to the System 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
aon exit, this contains the vector defining the first edge of the periodic box
bon exit, this contains the vector defining the second edge of the periodic box
con exit, this contains the vector defining the third edge of the periodic box

References simtk.openmm.openmm.stripUnits().

def getForce (   self,
  args 
)

getForce(System self, int index) -> Force getForce(System self, int index) -> Force

Get a writable reference to one of the Forces in this System.

Parameters
indexthe index of the Force to get

References simtk.openmm.openmm.stripUnits().

Referenced by System.getForces().

def getForces (   self)

Get the list of Forces in this System.

References System.getForce(), and System.getNumForces().

def getNumConstraints (   self)

getNumConstraints(System self) -> int

Get the number of distance constraints in this System.

References simtk.openmm.openmm.stripUnits().

def getNumForces (   self)

getNumForces(System self) -> int

Get the number of Force objects that have been added to the System.

References simtk.openmm.openmm.stripUnits().

Referenced by System.getForces().

def getNumParticles (   self)

getNumParticles(System self) -> int

Get the number of particles in this System.

References simtk.openmm.openmm.stripUnits().

def getParticleMass (   self,
  args 
)

getParticleMass(System self, int index) -> double

Get the mass (in atomic mass units) of a particle. If the mass is 0, Integrators will ignore the particle and not modify its position or velocity. This is most often used for virtual sites, but can also be used as a way to prevent a particle from moving.

Parameters
indexthe index of the particle for which to get the mass

References simtk.openmm.openmm.stripUnits().

def getVirtualSite (   self,
  args 
)

getVirtualSite(System self, int index) -> VirtualSite

Get VirtualSite object for a particle. If the particle is not a virtual site, this throws an exception.

Parameters
indexthe index of the particle to get

References simtk.openmm.openmm.stripUnits().

def isVirtualSite (   self,
  args 
)

isVirtualSite(System self, int index) -> bool

Get whether a particle is a VirtualSite.

Parameters
indexthe index of the particle to check

References simtk.openmm.openmm.stripUnits().

def setConstraintParameters (   self,
  args 
)

setConstraintParameters(System self, int index, int particle1, int particle2, double distance)

Set the parameters defining a distance constraint. Particles whose mass is 0 cannot participate in constraints.

Parameters
indexthe index of the constraint for which to set parameters
particle1the index of the first particle involved in the constraint
particle2the index of the second particle involved in the constraint
distancethe required distance between the two particles, measured in nm

References simtk.openmm.openmm.stripUnits().

def setDefaultPeriodicBoxVectors (   self,
  args 
)

setDefaultPeriodicBoxVectors(System self, Vec3 const & a, Vec3 const & b, Vec3 const & c)

Set the default values of the vectors defining the axes of the periodic box (measured in nm). Any newly created Context will have its box vectors set to these. They will affect any Force added to the System 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

References simtk.openmm.openmm.stripUnits().

def setParticleMass (   self,
  args 
)

setParticleMass(System self, int index, double mass)

Set the mass (in atomic mass units) of a particle. If the mass is 0, Integrators will ignore the particle and not modify its position or velocity. This is most often used for virtual sites, but can also be used as a way to prevent a particle from moving.

Parameters
indexthe index of the particle for which to set the mass
massthe mass of the particle

References simtk.openmm.openmm.stripUnits().

def setVirtualSite (   self,
  args 
)

setVirtualSite(System self, int index, VirtualSite virtualSite)

Set a particle to be a virtual site. The VirtualSite object should have been created on the heap with the "new" operator. The System takes over ownership of it, and deletes it when the System itself is deleted.

Parameters
indexthe index of the particle that should be treated as a virtual site
virtualSitea pointer to the VirtualSite object describing it

Member Data Documentation

this

Referenced by System.__init__().


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