This class implements an implicit solvation force using the GB/VI model. More...
Public Member Functions | |
def | getNumParticles |
getNumParticles(self) -> int | |
def | addParticle |
addParticle(self, double charge, double radius, double gamma) -> int | |
def | getParticleParameters |
getParticleParameters(self, int index) | |
def | setParticleParameters |
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) -> NonbondedMethod | |
def | setNonbondedMethod |
setNonbondedMethod(self, NonbondedMethod method) | |
def | getCutoffDistance |
getCutoffDistance(self) -> double | |
def | setCutoffDistance |
setCutoffDistance(self, double distance) | |
def | getBornRadiusScalingMethod |
getBornRadiusScalingMethod(self) -> BornRadiusScalingMethod | |
def | setBornRadiusScalingMethod |
setBornRadiusScalingMethod(self, BornRadiusScalingMethod 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) -> GBVIForce __init__(self, GBVIForce other) -> GBVIForce | |
def | __del__ |
__del__(self) | |
Public Attributes | |
this | |
Static Public Attributes | |
NoCutoff = _openmm.GBVIForce_NoCutoff | |
CutoffNonPeriodic = _openmm.GBVIForce_CutoffNonPeriodic | |
CutoffPeriodic = _openmm.GBVIForce_CutoffPeriodic | |
NoScaling = _openmm.GBVIForce_NoScaling | |
QuinticSpline = _openmm.GBVIForce_QuinticSpline |
This class implements an implicit solvation force using the GB/VI model.
To use this class, create a GBVIForce 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 a Context. After a particle has been added, you can modify its force field parameters by calling setParticleParameters().
def __del__ | ( | self | ) |
__del__(self)
Reimplemented from Force.
def __init__ | ( | self, | ||
args | ||||
) |
def addBond | ( | self, | ||
args | ||||
) |
addBond(self, int particle1, int particle2, double distance) -> int
Add a bond
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) -> int
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.
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 |
def getBondParameters | ( | self, | ||
args | ||||
) |
getBondParameters(self, int index)
Get the parameters defining a bond
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) -> BornRadiusScalingMethod
Get Born radius scaling method
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) -> NonbondedMethod
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.
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 |
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, measured in nm (~5.0)
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
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, measured in nm |
def setBornRadiusScalingMethod | ( | self, | ||
args | ||||
) |
setBornRadiusScalingMethod(self, BornRadiusScalingMethod method)
Set Born radius scaling method
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.
distance | the cutoff distance, measured in nm |
def setNonbondedMethod | ( | self, | ||
args | ||||
) |
setNonbondedMethod(self, NonbondedMethod 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)
Set the force field parameters for a particle.
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 |
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, measured in nm (~5.0)
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.
CutoffNonPeriodic = _openmm.GBVIForce_CutoffNonPeriodic [static] |
CutoffPeriodic = _openmm.GBVIForce_CutoffPeriodic [static] |
NoCutoff = _openmm.GBVIForce_NoCutoff [static] |
NoScaling = _openmm.GBVIForce_NoScaling [static] |
QuinticSpline = _openmm.GBVIForce_QuinticSpline [static] |