GBVISoftcoreForce Class Reference

This class implements an implicit solvation force using the GB/VI model. More...

Inheritance diagram for GBVISoftcoreForce:
Force

List of all members.

Public Member Functions

def getNumParticles
 getNumParticles(self) -> int
def addParticle
 addParticle(self, double charge, double radius, double gamma, double bornRadiusScaleFactor = 1.0) addParticle(self, double charge, double radius, double gamma)
def getParticleParameters
 getParticleParameters(self, int index)
def setParticleParameters
 setParticleParameters(self, int index, double charge, double radius, double gamma, double bornRadiusScaleFactor = 1.0) setParticleParameters(self, int index, double charge, double radius, double gamma)
def addBond
 addBond(self, int particle1, int particle2, double distance) -> int
def getBondParameters
 getBondParameters(self, int index)
def setBondParameters
 setBondParameters(self, int index, int particle1, int particle2, double bondLength)
def getNumBonds
 getNumBonds(self) -> int
def getSolventDielectric
 getSolventDielectric(self) -> double
def setSolventDielectric
 setSolventDielectric(self, double dielectric)
def getSoluteDielectric
 getSoluteDielectric(self) -> double
def setSoluteDielectric
 setSoluteDielectric(self, double dielectric)
def getNonbondedMethod
 getNonbondedMethod(self) -> NonbondedSoftcoreMethod
def setNonbondedMethod
 setNonbondedMethod(self, NonbondedSoftcoreMethod method)
def getCutoffDistance
 getCutoffDistance(self) -> double
def setCutoffDistance
 setCutoffDistance(self, double distance)
def getBornRadiusScalingMethod
 getBornRadiusScalingMethod(self) -> BornRadiusScalingSoftcoreMethod
def setBornRadiusScalingMethod
 setBornRadiusScalingMethod(self, BornRadiusScalingSoftcoreMethod method)
def getQuinticLowerLimitFactor
 getQuinticLowerLimitFactor(self) -> double
def setQuinticLowerLimitFactor
 setQuinticLowerLimitFactor(self, double quinticLowerLimitFactor)
def getQuinticUpperBornRadiusLimit
 getQuinticUpperBornRadiusLimit(self) -> double
def setQuinticUpperBornRadiusLimit
 setQuinticUpperBornRadiusLimit(self, double quinticUpperBornRadiusLimit)
def __init__
 __init__(self) -> GBVISoftcoreForce __init__(self, GBVISoftcoreForce other) -> GBVISoftcoreForce
def __del__
 __del__(self)

Public Attributes

 this

Static Public Attributes

 NoCutoff = _openmm.GBVISoftcoreForce_NoCutoff
 CutoffNonPeriodic = _openmm.GBVISoftcoreForce_CutoffNonPeriodic
 CutoffPeriodic = _openmm.GBVISoftcoreForce_CutoffPeriodic
 NoScaling = _openmm.GBVISoftcoreForce_NoScaling
 QuinticSpline = _openmm.GBVISoftcoreForce_QuinticSpline

Detailed Description

This class implements an implicit solvation force using the GB/VI model.

To use this class, create a GBVISoftcoreForce object, then call addParticle() once for each particle in the System to define its parameters. The number of particles for which you define GB/VI parameters must be exactly equal to the number of particles in the System, or else an exception will be thrown when you try to create an Context. After a particle has been added, you can modify its force field parameters by calling setParticleParameters().

If the System also contains a NonbondedForce, this force will use the cutoffs and periodic boundary conditions specified in it.


Member Function Documentation

def __del__ (   self  ) 

__del__(self)

Reimplemented from Force.

def __init__ (   self,
  args 
)

__init__(self) -> GBVISoftcoreForce __init__(self, GBVISoftcoreForce other) -> GBVISoftcoreForce

def addBond (   self,
  args 
)

addBond(self, int particle1, int particle2, double distance) -> int

Add a bond

Parameters:
particle1 the index of the first particle
particle2 the index of the second particle
distance the distance between the two particles, measured in nm
def addParticle (   self,
  args 
)

addParticle(self, double charge, double radius, double gamma, double bornRadiusScaleFactor = 1.0) addParticle(self, double charge, double radius, double gamma)

Add the GB/VI parameters for a particle. This should be called once for each particle in the System. When it is called for the i'th time, it specifies the parameters for the i'th particle.

Parameters:
charge the charge of the particle, measured in units of the proton charge
radius the GB/VI radius of the particle, measured in nm
gamma the gamma parameter
bornRadiusScaleFactor the Born radius scale factor (used for free energy calculations)
def getBondParameters (   self,
  args 
)

getBondParameters(self, int index)

Get the parameters defining a bond

Parameters:
index the index of the bond for which to get parameters
particle1 the index of the first particle involved in the bond
particle2 the index of the second particle involved in the bond
distance the distance between the two particles, measured in nm
def getBornRadiusScalingMethod (   self  ) 

getBornRadiusScalingMethod(self) -> BornRadiusScalingSoftcoreMethod

Get the method used for scaling Born radii.

def getCutoffDistance (   self  ) 

getCutoffDistance(self) -> double

Get the cutoff distance (in nm) being used for nonbonded interactions. If the NonbondedMethod in use is NoCutoff, this value will have no effect.

def getNonbondedMethod (   self  ) 

getNonbondedMethod(self) -> NonbondedSoftcoreMethod

Get the method used for handling long range nonbonded interactions.

def getNumBonds (   self  ) 

getNumBonds(self) -> int

Get number of bonds

def getNumParticles (   self  ) 

getNumParticles(self) -> int

Get the number of particles in the system.

def getParticleParameters (   self,
  args 
)

getParticleParameters(self, int index)

Get the force field parameters for a particle.

Parameters:
index the index of the particle for which to get parameters
charge the charge of the particle, measured in units of the proton charge
radius the GBSA radius of the particle, measured in nm
gamma the gamma parameter
bornRadiusScaleFactor the Born radius scale factor (used for free energy calculations)
def getQuinticLowerLimitFactor (   self  ) 

getQuinticLowerLimitFactor(self) -> double

Get the lower limit factor used in the quintic spline scaling method (typically 0.5-0.8)

def getQuinticUpperBornRadiusLimit (   self  ) 

getQuinticUpperBornRadiusLimit(self) -> double

Get the upper limit used in the quintic spline scaling method (typically 0.5-0.8)

def getSoluteDielectric (   self  ) 

getSoluteDielectric(self) -> double

Get the dielectric constant for the solute.

def getSolventDielectric (   self  ) 

getSolventDielectric(self) -> double

Get the dielectric constant for the solvent.

def setBondParameters (   self,
  args 
)

setBondParameters(self, int index, int particle1, int particle2, double bondLength)

Set 1-2 bonds

Parameters:
index index of the bond for which to set parameters
particle1 index of first atom in bond
particle2 index of second atom in bond
bondLength bond length
def setBornRadiusScalingMethod (   self,
  args 
)

setBornRadiusScalingMethod(self, BornRadiusScalingSoftcoreMethod method)

Set the method used for scaling Born radii.

def setCutoffDistance (   self,
  args 
)

setCutoffDistance(self, double distance)

Set the cutoff distance (in nm) being used for nonbonded interactions. If the NonbondedMethod in use is NoCutoff, this value will have no effect.

def setNonbondedMethod (   self,
  args 
)

setNonbondedMethod(self, NonbondedSoftcoreMethod method)

Set the method used for handling long range nonbonded interactions.

def setParticleParameters (   self,
  args 
)

setParticleParameters(self, int index, double charge, double radius, double gamma, double bornRadiusScaleFactor = 1.0) setParticleParameters(self, int index, double charge, double radius, double gamma)

Set the force field parameters for a particle.

Parameters:
index the index of the particle for which to set parameters
charge the charge of the particle, measured in units of the proton charge
radius the GB/VI radius of the particle, measured in nm
gamma the gamma parameter
bornRadiusScaleFactor the Born radius scale factor (used for free energy calculations)
def setQuinticLowerLimitFactor (   self,
  args 
)

setQuinticLowerLimitFactor(self, double quinticLowerLimitFactor)

Set the lower limit factor used in the quintic spline scaling method (typically 0.5-0.8)

def setQuinticUpperBornRadiusLimit (   self,
  args 
)

setQuinticUpperBornRadiusLimit(self, double quinticUpperBornRadiusLimit)

Set the upper limit used in the quintic spline scaling method (typically 0.008)

def setSoluteDielectric (   self,
  args 
)

setSoluteDielectric(self, double dielectric)

Set the dielectric constant for the solute.

def setSolventDielectric (   self,
  args 
)

setSolventDielectric(self, double dielectric)

Set the dielectric constant for the solvent.


Member Data Documentation

CutoffNonPeriodic = _openmm.GBVISoftcoreForce_CutoffNonPeriodic [static]
CutoffPeriodic = _openmm.GBVISoftcoreForce_CutoffPeriodic [static]
NoCutoff = _openmm.GBVISoftcoreForce_NoCutoff [static]
NoScaling = _openmm.GBVISoftcoreForce_NoScaling [static]
QuinticSpline = _openmm.GBVISoftcoreForce_QuinticSpline [static]

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

Generated by  doxygen 1.6.2