OpenMM
|
This class implements an implicit solvation force using the GB/VI model. More...
Public Member Functions | |
def | getNumParticles |
getNumParticles(GBVIForce self) -> int | |
def | addParticle |
addParticle(GBVIForce self, double charge, double radius, double gamma) -> int | |
def | getParticleParameters |
getParticleParameters(GBVIForce self, int index) | |
def | setParticleParameters |
setParticleParameters(GBVIForce self, int index, double charge, double radius, double gamma) | |
def | addBond |
addBond(GBVIForce self, int particle1, int particle2, double distance) -> int | |
def | getBondParameters |
getBondParameters(GBVIForce self, int index) | |
def | setBondParameters |
setBondParameters(GBVIForce self, int index, int particle1, int particle2, double bondLength) | |
def | getNumBonds |
getNumBonds(GBVIForce self) -> int | |
def | getSolventDielectric |
getSolventDielectric(GBVIForce self) -> double | |
def | setSolventDielectric |
setSolventDielectric(GBVIForce self, double dielectric) | |
def | getSoluteDielectric |
getSoluteDielectric(GBVIForce self) -> double | |
def | setSoluteDielectric |
setSoluteDielectric(GBVIForce self, double dielectric) | |
def | getNonbondedMethod |
getNonbondedMethod(GBVIForce self) -> OpenMM::GBVIForce::NonbondedMethod | |
def | setNonbondedMethod |
setNonbondedMethod(GBVIForce self, OpenMM::GBVIForce::NonbondedMethod method) | |
def | getCutoffDistance |
getCutoffDistance(GBVIForce self) -> double | |
def | setCutoffDistance |
setCutoffDistance(GBVIForce self, double distance) | |
def | getBornRadiusScalingMethod |
getBornRadiusScalingMethod(GBVIForce self) -> OpenMM::GBVIForce::BornRadiusScalingMethod | |
def | setBornRadiusScalingMethod |
setBornRadiusScalingMethod(GBVIForce self, OpenMM::GBVIForce::BornRadiusScalingMethod method) | |
def | getQuinticLowerLimitFactor |
getQuinticLowerLimitFactor(GBVIForce self) -> double | |
def | setQuinticLowerLimitFactor |
setQuinticLowerLimitFactor(GBVIForce self, double quinticLowerLimitFactor) | |
def | getQuinticUpperBornRadiusLimit |
getQuinticUpperBornRadiusLimit(GBVIForce self) -> double | |
def | setQuinticUpperBornRadiusLimit |
setQuinticUpperBornRadiusLimit(GBVIForce self, double quinticUpperBornRadiusLimit) | |
def | __init__ |
init(OpenMM::GBVIForce self) -> GBVIForce init(OpenMM::GBVIForce self, GBVIForce other) -> GBVIForce | |
def | __del__ |
del(OpenMM::GBVIForce self) | |
Public Member Functions inherited from Force | |
def | __init__ |
def | __del__ |
del(OpenMM::Force self) | |
def | getForceGroup |
getForceGroup(Force self) -> int | |
def | setForceGroup |
setForceGroup(Force self, int group) | |
def | __copy__ |
def | __deepcopy__ |
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 __init__ | ( | self, | |
args | |||
) |
init(OpenMM::GBVIForce self) -> GBVIForce init(OpenMM::GBVIForce self, GBVIForce other) -> GBVIForce
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(OpenMM::GBVIForce self)
def addBond | ( | self, | |
args | |||
) |
addBond(GBVIForce 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(GBVIForce 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(GBVIForce 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(GBVIForce self) -> OpenMM::GBVIForce::BornRadiusScalingMethod
Get Born radius scaling method
def getCutoffDistance | ( | self | ) |
getCutoffDistance(GBVIForce 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(GBVIForce self) -> OpenMM::GBVIForce::NonbondedMethod
Get the method used for handling long range nonbonded interactions.
def getNumBonds | ( | self | ) |
getNumBonds(GBVIForce self) -> int
Get number of bonds
def getNumParticles | ( | self | ) |
getNumParticles(GBVIForce self) -> int
Get the number of particles in the system.
def getParticleParameters | ( | self, | |
args | |||
) |
getParticleParameters(GBVIForce 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(GBVIForce self) -> double
Get the lower limit factor used in the quintic spline scaling method (typically 0.5-0.8)
def getQuinticUpperBornRadiusLimit | ( | self | ) |
getQuinticUpperBornRadiusLimit(GBVIForce self) -> double
Get the upper limit used in the quintic spline scaling method, measured in nm (~5.0)
def getSoluteDielectric | ( | self | ) |
getSoluteDielectric(GBVIForce self) -> double
Get the dielectric constant for the solute.
def getSolventDielectric | ( | self | ) |
getSolventDielectric(GBVIForce self) -> double
Get the dielectric constant for the solvent.
def setBondParameters | ( | self, | |
args | |||
) |
setBondParameters(GBVIForce 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(GBVIForce self, OpenMM::GBVIForce::BornRadiusScalingMethod method)
Set Born radius scaling method
def setCutoffDistance | ( | self, | |
args | |||
) |
setCutoffDistance(GBVIForce 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(GBVIForce self, OpenMM::GBVIForce::NonbondedMethod method)
Set the method used for handling long range nonbonded interactions.
def setParticleParameters | ( | self, | |
args | |||
) |
setParticleParameters(GBVIForce 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(GBVIForce 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(GBVIForce self, double quinticUpperBornRadiusLimit)
Set the upper limit used in the quintic spline scaling method, measured in nm (~5.0)
def setSoluteDielectric | ( | self, | |
args | |||
) |
setSoluteDielectric(GBVIForce self, double dielectric)
Set the dielectric constant for the solute.
def setSolventDielectric | ( | self, | |
args | |||
) |
setSolventDielectric(GBVIForce self, double dielectric)
Set the dielectric constant for the solvent.
|
static |
|
static |
|
static |
|
static |
|
static |
this |