OpenMM
 All Classes Namespaces Functions Variables Pages
MonteCarloMembraneBarostat Class Reference

This is a Monte Carlo barostat designed specifically for membrane simulations. More...

+ Inheritance diagram for MonteCarloMembraneBarostat:

Public Member Functions

def Pressure
 Pressure() -> std::string const &. More...
 
def SurfaceTension
 SurfaceTension() -> std::string const &. More...
 
def getDefaultPressure
 getDefaultPressure(MonteCarloMembraneBarostat self) -> double More...
 
def setDefaultPressure
 setDefaultPressure(MonteCarloMembraneBarostat self, double pressure) More...
 
def getDefaultSurfaceTension
 getDefaultSurfaceTension(MonteCarloMembraneBarostat self) -> double More...
 
def setDefaultSurfaceTension
 setDefaultSurfaceTension(MonteCarloMembraneBarostat self, double surfaceTension) More...
 
def getFrequency
 getFrequency(MonteCarloMembraneBarostat self) -> int More...
 
def setFrequency
 setFrequency(MonteCarloMembraneBarostat self, int freq) More...
 
def getTemperature
 getTemperature(MonteCarloMembraneBarostat self) -> double More...
 
def setTemperature
 setTemperature(MonteCarloMembraneBarostat self, double temp) More...
 
def getXYMode
 getXYMode(MonteCarloMembraneBarostat self) -> OpenMM::MonteCarloMembraneBarostat::XYMode More...
 
def setXYMode
 setXYMode(MonteCarloMembraneBarostat self, OpenMM::MonteCarloMembraneBarostat::XYMode mode) More...
 
def getZMode
 getZMode(MonteCarloMembraneBarostat self) -> OpenMM::MonteCarloMembraneBarostat::ZMode More...
 
def setZMode
 setZMode(MonteCarloMembraneBarostat self, OpenMM::MonteCarloMembraneBarostat::ZMode mode) More...
 
def getRandomNumberSeed
 getRandomNumberSeed(MonteCarloMembraneBarostat self) -> int More...
 
def setRandomNumberSeed
 setRandomNumberSeed(MonteCarloMembraneBarostat self, int seed) More...
 
def usesPeriodicBoundaryConditions
 usesPeriodicBoundaryConditions(MonteCarloMembraneBarostat self) -> bool More...
 
def __init__
 init(OpenMM::MonteCarloMembraneBarostat self, double defaultPressure, double defaultSurfaceTension, double temperature, OpenMM::MonteCarloMembraneBarostat::XYMode xymode, OpenMM::MonteCarloMembraneBarostat::ZMode zmode, int frequency=25) -> MonteCarloMembraneBarostat init(OpenMM::MonteCarloMembraneBarostat self, MonteCarloMembraneBarostat other) -> MonteCarloMembraneBarostat More...
 
def __del__
 del(OpenMM::MonteCarloMembraneBarostat self) More...
 
- Public Member Functions inherited from Force
def __init__
 
def __del__
 del(OpenMM::Force self) More...
 
def getForceGroup
 getForceGroup(Force self) -> int More...
 
def setForceGroup
 setForceGroup(Force self, int group) More...
 
def usesPeriodicBoundaryConditions
 usesPeriodicBoundaryConditions(Force self) -> bool More...
 
def __copy__
 
def __deepcopy__
 

Public Attributes

 this
 

Static Public Attributes

 XYIsotropic = _openmm.MonteCarloMembraneBarostat_XYIsotropic
 
 XYAnisotropic = _openmm.MonteCarloMembraneBarostat_XYAnisotropic
 
 ZFree = _openmm.MonteCarloMembraneBarostat_ZFree
 
 ZFixed = _openmm.MonteCarloMembraneBarostat_ZFixed
 
 ConstantVolume = _openmm.MonteCarloMembraneBarostat_ConstantVolume
 

Detailed Description

This is a Monte Carlo barostat designed specifically for membrane simulations.

It assumes the membrane lies in the XY plane. The Monte Carlo acceptance criterion includes a term to model isotropic pressure, which depends on the volume of the periodic box, and a second term to model surface tension, which depends on the cross sectional area of the box in the XY plane. Note that pressure and surface tension are defined with opposite senses: a larger pressure tends to make the box smaller, but a larger surface tension tends to make the box larger.

There are options for configuring exactly how the various box dimensions are allowed to change:

  • The X and Y axes may be treated isotropically, in which case they always scale by the same amount and remain in proportion to each other; or they may be treated anisotropically, in which case they can vary independently of each other.

  • The Z axis can be allowed to vary independently of the other axes; or held fixed; or constrained to vary in inverse proportion to the other two axes, so that the total box volume remains fixed.

This class assumes the simulation is also being run at constant temperature, and requires you to specify the system temperature (since it affects the acceptance probability for Monte Carlo moves). It does not actually perform temperature regulation, however. You must use another mechanism along with it to maintain the temperature, such as LangevinIntegrator or AndersenThermostat.

Constructor & Destructor Documentation

def __init__ (   self,
  args 
)

init(OpenMM::MonteCarloMembraneBarostat self, double defaultPressure, double defaultSurfaceTension, double temperature, OpenMM::MonteCarloMembraneBarostat::XYMode xymode, OpenMM::MonteCarloMembraneBarostat::ZMode zmode, int frequency=25) -> MonteCarloMembraneBarostat init(OpenMM::MonteCarloMembraneBarostat self, MonteCarloMembraneBarostat other) -> MonteCarloMembraneBarostat

Create a MonteCarloMembraneBarostat.

Parameters
defaultPressurethe default pressure acting on the system (in bar)
defaultSurfaceTensionthe default surface tension acting on the system (in bar*nm)
temperaturethe temperature at which the system is being maintained (in Kelvin)
xymodethe mode specifying the behavior of the X and Y axes
zmodethe mode specifying the behavior of the Z axis
frequencythe frequency at which Monte Carlo volume changes should be attempted (in time steps)

References simtk.openmm.openmm.stripUnits().

def __del__ (   self)

del(OpenMM::MonteCarloMembraneBarostat self)

References simtk.openmm.openmm.stripUnits().

Member Function Documentation

def getDefaultPressure (   self,
  args 
)

getDefaultPressure(MonteCarloMembraneBarostat self) -> double

Get the default pressure acting on the system (in bar).

References simtk.openmm.openmm.stripUnits().

def getDefaultSurfaceTension (   self,
  args 
)

getDefaultSurfaceTension(MonteCarloMembraneBarostat self) -> double

Get the default surface tension acting on the system (in bar*nm).

References simtk.openmm.openmm.stripUnits().

def getFrequency (   self,
  args 
)

getFrequency(MonteCarloMembraneBarostat self) -> int

Get the frequency (in time steps) at which Monte Carlo volume changes should be attempted. If this is set to 0, the barostat is disabled.

References simtk.openmm.openmm.stripUnits().

def getRandomNumberSeed (   self,
  args 
)

getRandomNumberSeed(MonteCarloMembraneBarostat self) -> int

Get the random number seed. See setRandomNumberSeed() for details.

References simtk.openmm.openmm.stripUnits().

def getTemperature (   self,
  args 
)

getTemperature(MonteCarloMembraneBarostat self) -> double

Get the temperature at which the system is being maintained, measured in Kelvin.

References simtk.openmm.openmm.stripUnits().

def getXYMode (   self,
  args 
)

getXYMode(MonteCarloMembraneBarostat self) -> OpenMM::MonteCarloMembraneBarostat::XYMode

Get the mode specifying the behavior of the X and Y axes.

References simtk.openmm.openmm.stripUnits().

def getZMode (   self,
  args 
)

getZMode(MonteCarloMembraneBarostat self) -> OpenMM::MonteCarloMembraneBarostat::ZMode

Get the mode specifying the behavior of the Z axis.

References simtk.openmm.openmm.stripUnits().

def Pressure (   args)

Pressure() -> std::string const &.

This is the name of the parameter which stores the current pressure acting on the system (in bar).

References simtk.openmm.openmm.stripUnits().

def setDefaultPressure (   self,
  args 
)

setDefaultPressure(MonteCarloMembraneBarostat self, double pressure)

Set the default pressure acting on the system. This will affect any new Contexts you create, but not ones that already exist.

Parameters
pressurethe default pressure acting on the system, measured in bar.

References simtk.openmm.openmm.stripUnits().

def setDefaultSurfaceTension (   self,
  args 
)

setDefaultSurfaceTension(MonteCarloMembraneBarostat self, double surfaceTension)

Set the default surface tension acting on the system. This will affect any new Contexts you create, but not ones that already exist.

Parameters
surfaceTensionthe default surface tension acting on the system, measured in bar.

References simtk.openmm.openmm.stripUnits().

def setFrequency (   self,
  args 
)

setFrequency(MonteCarloMembraneBarostat self, int freq)

Set the frequency (in time steps) at which Monte Carlo volume changes should be attempted. If this is set to 0, the barostat is disabled.

References simtk.openmm.openmm.stripUnits().

def setRandomNumberSeed (   self,
  args 
)

setRandomNumberSeed(MonteCarloMembraneBarostat self, int seed)

Set the random number seed. It is guaranteed that if two simulations are run with different random number seeds, the sequence of Monte Carlo steps will be different. On the other hand, no guarantees are made about the behavior of simulations that use the same seed. In particular, Platforms are permitted to use non-deterministic algorithms which produce different results on successive runs, even if those runs were initialized identically.

If seed is set to 0 (which is the default value assigned), a unique seed is chosen when a Context is created from this Force. This is done to ensure that each Context receives unique random seeds without you needing to set them explicitly.

References simtk.openmm.openmm.stripUnits().

def setTemperature (   self,
  args 
)

setTemperature(MonteCarloMembraneBarostat self, double temp)

Set the temperature at which the system is being maintained.

Parameters
tempthe system temperature, measured in Kelvin.

References simtk.openmm.openmm.stripUnits().

def setXYMode (   self,
  args 
)

setXYMode(MonteCarloMembraneBarostat self, OpenMM::MonteCarloMembraneBarostat::XYMode mode)

Set the mode specifying the behavior of the X and Y axes.

References simtk.openmm.openmm.stripUnits().

def setZMode (   self,
  args 
)

setZMode(MonteCarloMembraneBarostat self, OpenMM::MonteCarloMembraneBarostat::ZMode mode)

Set the mode specifying the behavior of the Z axis.

References simtk.openmm.openmm.stripUnits().

def SurfaceTension (   args)

SurfaceTension() -> std::string const &.

This is the name of the parameter which stores the current surface tension acting on the system (in bar*nm).

References simtk.openmm.openmm.stripUnits().

def usesPeriodicBoundaryConditions (   self,
  args 
)

usesPeriodicBoundaryConditions(MonteCarloMembraneBarostat self) -> bool

Returns whether or not this force makes use of periodic boundary conditions.

References simtk.openmm.openmm.stripUnits().

Member Data Documentation

ConstantVolume = _openmm.MonteCarloMembraneBarostat_ConstantVolume
static
this

Referenced by System.__init__().

XYAnisotropic = _openmm.MonteCarloMembraneBarostat_XYAnisotropic
static
XYIsotropic = _openmm.MonteCarloMembraneBarostat_XYIsotropic
static
ZFixed = _openmm.MonteCarloMembraneBarostat_ZFixed
static
ZFree = _openmm.MonteCarloMembraneBarostat_ZFree
static

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