OpenMM
 All Classes Namespaces Functions Variables Pages
GBSAOBCForce Class Reference

This class implements an implicit solvation force using the GBSA-OBC model. More...

+ Inheritance diagram for GBSAOBCForce:

Public Member Functions

def getNumParticles
 getNumParticles(GBSAOBCForce self) -> int More...
 
def addParticle
 addParticle(GBSAOBCForce self, double charge, double radius, double scalingFactor) -> int More...
 
def getParticleParameters
 getParticleParameters(GBSAOBCForce self, int index) More...
 
def setParticleParameters
 setParticleParameters(GBSAOBCForce self, int index, double charge, double radius, double scalingFactor) More...
 
def getSolventDielectric
 getSolventDielectric(GBSAOBCForce self) -> double More...
 
def setSolventDielectric
 setSolventDielectric(GBSAOBCForce self, double dielectric) More...
 
def getSoluteDielectric
 getSoluteDielectric(GBSAOBCForce self) -> double More...
 
def setSoluteDielectric
 setSoluteDielectric(GBSAOBCForce self, double dielectric) More...
 
def getSurfaceAreaEnergy
 getSurfaceAreaEnergy(GBSAOBCForce self) -> double More...
 
def setSurfaceAreaEnergy
 setSurfaceAreaEnergy(GBSAOBCForce self, double energy) More...
 
def getNonbondedMethod
 getNonbondedMethod(GBSAOBCForce self) -> OpenMM::GBSAOBCForce::NonbondedMethod More...
 
def setNonbondedMethod
 setNonbondedMethod(GBSAOBCForce self, OpenMM::GBSAOBCForce::NonbondedMethod method) More...
 
def getCutoffDistance
 getCutoffDistance(GBSAOBCForce self) -> double More...
 
def setCutoffDistance
 setCutoffDistance(GBSAOBCForce self, double distance) More...
 
def updateParametersInContext
 updateParametersInContext(GBSAOBCForce self, Context context) More...
 
def usesPeriodicBoundaryConditions
 usesPeriodicBoundaryConditions(GBSAOBCForce self) -> bool More...
 
def __init__
 init(OpenMM::GBSAOBCForce self) -> GBSAOBCForce init(OpenMM::GBSAOBCForce self, GBSAOBCForce other) -> GBSAOBCForce More...
 
def __del__
 del(OpenMM::GBSAOBCForce 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

 NoCutoff = _openmm.GBSAOBCForce_NoCutoff
 
 CutoffNonPeriodic = _openmm.GBSAOBCForce_CutoffNonPeriodic
 
 CutoffPeriodic = _openmm.GBSAOBCForce_CutoffPeriodic
 

Detailed Description

This class implements an implicit solvation force using the GBSA-OBC model.

To use this class, create a GBSAOBCForce object, then call addParticle() once for each particle in the System to define its parameters. The number of particles for which you define GBSA 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 a Context. After a particle has been added, you can modify its force field parameters by calling setParticleParameters(). This will have no effect on Contexts that already exist unless you call updateParametersInContext().

When using this Force, the System should also include a NonbondedForce, and both objects must specify identical charges for all particles. Otherwise, the results will not be correct. Furthermore, if the nonbonded method is set to CutoffNonPeriodic or CutoffPeriodic, you should call setReactionFieldDielectric(1.0) on the NonbondedForce to turn off the reaction field approximation, which does not produce correct results when combined with GBSA.

Constructor & Destructor Documentation

def __init__ (   self,
  args 
)

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

Create a GBSAOBCForce.

References simtk.openmm.openmm.stripUnits().

def __del__ (   self)

del(OpenMM::GBSAOBCForce self)

References simtk.openmm.openmm.stripUnits().

Member Function Documentation

def addParticle (   self,
  args 
)

addParticle(GBSAOBCForce self, double charge, double radius, double scalingFactor) -> int

Add the GBSA 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
chargethe charge of the particle, measured in units of the proton charge
radiusthe GBSA radius of the particle, measured in nm
scalingFactorthe OBC scaling factor for the particle

References simtk.openmm.openmm.stripUnits().

Referenced by NonbondedForce.addParticle_usingRVdw().

def getCutoffDistance (   self,
  args 
)

getCutoffDistance(GBSAOBCForce 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.

References simtk.openmm.openmm.stripUnits().

def getNonbondedMethod (   self,
  args 
)

getNonbondedMethod(GBSAOBCForce self) -> OpenMM::GBSAOBCForce::NonbondedMethod

Get the method used for handling long range nonbonded interactions.

References simtk.openmm.openmm.stripUnits().

def getNumParticles (   self,
  args 
)

getNumParticles(GBSAOBCForce self) -> int

Get the number of particles in the system.

References simtk.openmm.openmm.stripUnits().

def getParticleParameters (   self,
  args 
)

getParticleParameters(GBSAOBCForce self, int index)

Get the force field parameters for a particle.

Parameters
indexthe index of the particle for which to get parameters
chargethe charge of the particle, measured in units of the proton charge
radiusthe GBSA radius of the particle, measured in nm
scalingFactorthe OBC scaling factor for the particle

References simtk.openmm.openmm.stripUnits().

def getSoluteDielectric (   self,
  args 
)

getSoluteDielectric(GBSAOBCForce self) -> double

Get the dielectric constant for the solute.

References simtk.openmm.openmm.stripUnits().

def getSolventDielectric (   self,
  args 
)

getSolventDielectric(GBSAOBCForce self) -> double

Get the dielectric constant for the solvent.

References simtk.openmm.openmm.stripUnits().

def getSurfaceAreaEnergy (   self,
  args 
)

getSurfaceAreaEnergy(GBSAOBCForce self) -> double

Get the energy scale for the surface energy term, measured in kJ/mol/nm^2.

References simtk.openmm.openmm.stripUnits().

def setCutoffDistance (   self,
  args 
)

setCutoffDistance(GBSAOBCForce 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.

Parameters
distancethe cutoff distance, measured in nm

References simtk.openmm.openmm.stripUnits().

def setNonbondedMethod (   self,
  args 
)

setNonbondedMethod(GBSAOBCForce self, OpenMM::GBSAOBCForce::NonbondedMethod method)

Set the method used for handling long range nonbonded interactions.

References simtk.openmm.openmm.stripUnits().

def setParticleParameters (   self,
  args 
)

setParticleParameters(GBSAOBCForce self, int index, double charge, double radius, double scalingFactor)

Set the force field parameters for a particle.

Parameters
indexthe index of the particle for which to set parameters
chargethe charge of the particle, measured in units of the proton charge
radiusthe GBSA radius of the particle, measured in nm
scalingFactorthe OBC scaling factor for the particle

References simtk.openmm.openmm.stripUnits().

def setSoluteDielectric (   self,
  args 
)

setSoluteDielectric(GBSAOBCForce self, double dielectric)

Set the dielectric constant for the solute.

References simtk.openmm.openmm.stripUnits().

def setSolventDielectric (   self,
  args 
)

setSolventDielectric(GBSAOBCForce self, double dielectric)

Set the dielectric constant for the solvent.

References simtk.openmm.openmm.stripUnits().

def setSurfaceAreaEnergy (   self,
  args 
)

setSurfaceAreaEnergy(GBSAOBCForce self, double energy)

Set the energy scale for the surface energy term, measured in kJ/mol/nm^2.

References simtk.openmm.openmm.stripUnits().

def updateParametersInContext (   self,
  args 
)

updateParametersInContext(GBSAOBCForce self, Context context)

Update the particle parameters in a Context to match those stored in this Force object. This method provides an efficient method to update certain parameters in an existing Context without needing to reinitialize it. Simply call setParticleParameters() to modify this object's parameters, then call updateParametersInContext() to copy them over to the Context.

The only information this method updates is the values of per-particle parameters. All other aspects of the Force (the nonbonded method, the cutoff distance, etc.) are unaffected and can only be changed by reinitializing the Context. Furthermore, this method cannot be used to add new particles, only to change the parameters of existing ones.

References simtk.openmm.openmm.stripUnits().

def usesPeriodicBoundaryConditions (   self,
  args 
)

usesPeriodicBoundaryConditions(GBSAOBCForce self) -> bool

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

References simtk.openmm.openmm.stripUnits().

Member Data Documentation

CutoffNonPeriodic = _openmm.GBSAOBCForce_CutoffNonPeriodic
static
CutoffPeriodic = _openmm.GBSAOBCForce_CutoffPeriodic
static
NoCutoff = _openmm.GBSAOBCForce_NoCutoff
static
this

Referenced by System.__init__().


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