OpenMM
|
This class implements an implicit solvation force using the GBSA-OBC model. More...
Public Member Functions | |
def | getNumParticles |
getNumParticles(GBSAOBCForce self) -> int | |
def | addParticle |
addParticle(GBSAOBCForce self, double charge, double radius, double scalingFactor) -> int | |
def | getParticleParameters |
getParticleParameters(GBSAOBCForce self, int index) | |
def | setParticleParameters |
setParticleParameters(GBSAOBCForce self, int index, double charge, double radius, double scalingFactor) | |
def | getSolventDielectric |
getSolventDielectric(GBSAOBCForce self) -> double | |
def | setSolventDielectric |
setSolventDielectric(GBSAOBCForce self, double dielectric) | |
def | getSoluteDielectric |
getSoluteDielectric(GBSAOBCForce self) -> double | |
def | setSoluteDielectric |
setSoluteDielectric(GBSAOBCForce self, double dielectric) | |
def | getNonbondedMethod |
getNonbondedMethod(GBSAOBCForce self) -> OpenMM::GBSAOBCForce::NonbondedMethod | |
def | setNonbondedMethod |
setNonbondedMethod(GBSAOBCForce self, OpenMM::GBSAOBCForce::NonbondedMethod method) | |
def | getCutoffDistance |
getCutoffDistance(GBSAOBCForce self) -> double | |
def | setCutoffDistance |
setCutoffDistance(GBSAOBCForce self, double distance) | |
def | updateParametersInContext |
updateParametersInContext(GBSAOBCForce self, Context context) | |
def | __init__ |
init(OpenMM::GBSAOBCForce self) -> GBSAOBCForce init(OpenMM::GBSAOBCForce self, GBSAOBCForce other) -> GBSAOBCForce | |
def | __del__ |
del(OpenMM::GBSAOBCForce self) | |
![]() | |
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.GBSAOBCForce_NoCutoff | |
CutoffNonPeriodic _openmm.GBSAOBCForce_CutoffNonPeriodic | |
CutoffPeriodic _openmm.GBSAOBCForce_CutoffPeriodic | |
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().
def __init__ | ( | self, | |
args | |||
) |
init(OpenMM::GBSAOBCForce self) -> GBSAOBCForce init(OpenMM::GBSAOBCForce self, GBSAOBCForce other) -> GBSAOBCForce
def __del__ | ( | self | ) |
del(OpenMM::GBSAOBCForce self)
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.
charge | the charge of the particle, measured in units of the proton charge |
radius | the GBSA radius of the particle, measured in nm |
scalingFactor | the OBC scaling factor for the particle |
def getCutoffDistance | ( | self | ) |
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.
def getNonbondedMethod | ( | self | ) |
getNonbondedMethod(GBSAOBCForce self) -> OpenMM::GBSAOBCForce::NonbondedMethod
Get the method used for handling long range nonbonded interactions.
def getNumParticles | ( | self | ) |
getNumParticles(GBSAOBCForce self) -> int
Get the number of particles in the system.
def getParticleParameters | ( | self, | |
args | |||
) |
getParticleParameters(GBSAOBCForce 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 |
scalingFactor | the OBC scaling factor for the particle |
def getSoluteDielectric | ( | self | ) |
getSoluteDielectric(GBSAOBCForce self) -> double
Get the dielectric constant for the solute.
def getSolventDielectric | ( | self | ) |
getSolventDielectric(GBSAOBCForce self) -> double
Get the dielectric constant for the solvent.
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.
distance | the cutoff distance, measured in nm |
def setNonbondedMethod | ( | self, | |
args | |||
) |
setNonbondedMethod(GBSAOBCForce self, OpenMM::GBSAOBCForce::NonbondedMethod method)
Set the method used for handling long range nonbonded interactions.
def setParticleParameters | ( | self, | |
args | |||
) |
setParticleParameters(GBSAOBCForce self, int index, double charge, double radius, double scalingFactor)
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 GBSA radius of the particle, measured in nm |
scalingFactor | the OBC scaling factor for the particle |
def setSoluteDielectric | ( | self, | |
args | |||
) |
setSoluteDielectric(GBSAOBCForce self, double dielectric)
Set the dielectric constant for the solute.
def setSolventDielectric | ( | self, | |
args | |||
) |
setSolventDielectric(GBSAOBCForce self, double dielectric)
Set the dielectric constant for the solvent.
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 updateParametersInState() 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.
|
static |
|
static |
|
static |
this |