OpenMM
GBSAOBCForce Class Reference

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

+ Inheritance diagram for GBSAOBCForce:

List of all members.

Public Member Functions

def getNumParticles
 getNumParticles(self) -> int
def addParticle
 addParticle(self, charge, radius, scalingFactor) -> int
def getParticleParameters
 Get the force field parameters for a particle.
def setParticleParameters
 Set the force field parameters for a particle.
def getSolventDielectric
 getSolventDielectric(self) -> double
def setSolventDielectric
 Set the dielectric constant for the solvent.
def getSoluteDielectric
 getSoluteDielectric(self) -> double
def setSoluteDielectric
 Set the dielectric constant for the solute.
def getSurfaceAreaEnergy
 getSurfaceAreaEnergy(self) -> double
def setSurfaceAreaEnergy
 Set the energy scale for the surface energy term, measured in kJ/mol/nm^2.
def getNonbondedMethod
 getNonbondedMethod(self) -> OpenMM::GBSAOBCForce::NonbondedMethod
def setNonbondedMethod
 Set the method used for handling long range nonbonded interactions.
def getCutoffDistance
 getCutoffDistance(self) -> double
def setCutoffDistance
 Set the cutoff distance (in nm) being used for nonbonded interactions.
def updateParametersInContext
 Update the particle parameters in a Context to match those stored in this Force object.
def usesPeriodicBoundaryConditions
 usesPeriodicBoundaryConditions(self) -> bool
def __init__
 __init__(self) -> GBSAOBCForce __init__(self, other) -> GBSAOBCForce

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__(self) -> GBSAOBCForce __init__(self, other) -> GBSAOBCForce

Create a GBSAOBCForce.


Member Function Documentation

def addParticle (   self,
  charge,
  radius,
  scalingFactor 
)

addParticle(self, charge, radius, 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:
charge(double) the charge of the particle, measured in units of the proton charge
radius(double) the GBSA radius of the particle, measured in nm
scalingFactor(double) the OBC scaling factor for the particle
Returns:
(int) the index of the particle that was added
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.

Returns:
(double) the cutoff distance, measured in nm
def getNonbondedMethod (   self)

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

Get the method used for handling long range nonbonded interactions.

def getNumParticles (   self)

getNumParticles(self) -> int

Get the number of particles in the system.

def getParticleParameters (   self,
  index 
)

Get the force field parameters for a particle.

Parameters:
index(int) the index of the particle for which to get parameters
Returns:
(double) the charge of the particle, measured in units of the proton charge
(double) the GBSA radius of the particle, measured in nm
(double) the OBC scaling factor for the particle
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 getSurfaceAreaEnergy (   self)

getSurfaceAreaEnergy(self) -> double

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

def setCutoffDistance (   self,
  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:
distance(double) the cutoff distance, measured in nm
def setNonbondedMethod (   self,
  method 
)

Set the method used for handling long range nonbonded interactions.

def setParticleParameters (   self,
  index,
  charge,
  radius,
  scalingFactor 
)

Set the force field parameters for a particle.

Parameters:
index(int) the index of the particle for which to set parameters
charge(double) the charge of the particle, measured in units of the proton charge
radius(double) the GBSA radius of the particle, measured in nm
scalingFactor(double) the OBC scaling factor for the particle
def setSoluteDielectric (   self,
  dielectric 
)

Set the dielectric constant for the solute.

def setSolventDielectric (   self,
  dielectric 
)

Set the dielectric constant for the solvent.

def setSurfaceAreaEnergy (   self,
  energy 
)

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

def updateParametersInContext (   self,
  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.

usesPeriodicBoundaryConditions(self) -> bool

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

Returns:
(bool) true if force uses PBC and false otherwise

Reimplemented from Force.


Member Data Documentation

CutoffNonPeriodic = _openmm.GBSAOBCForce_CutoffNonPeriodic [static]
CutoffPeriodic = _openmm.GBSAOBCForce_CutoffPeriodic [static]
NoCutoff = _openmm.GBSAOBCForce_NoCutoff [static]

Reimplemented from Force.


The documentation for this class was generated from the following file:
 All Classes Functions Variables