OpenMM
|
This class implements a buffered 14-7 potential used to model van der Waals forces. More...
Public Member Functions | |
def | getNumParticles |
getNumParticles(AmoebaVdwForce self) -> int | |
def | setParticleParameters |
setParticleParameters(AmoebaVdwForce self, int particleIndex, int parentIndex, double sigma, double epsilon, double reductionFactor) | |
def | getParticleParameters |
getParticleParameters(AmoebaVdwForce self, int particleIndex) | |
def | addParticle |
addParticle(AmoebaVdwForce self, int parentIndex, double sigma, double epsilon, double reductionFactor) -> int | |
def | setSigmaCombiningRule |
setSigmaCombiningRule(AmoebaVdwForce self, std::string const & sigmaCombiningRule) | |
def | getSigmaCombiningRule |
getSigmaCombiningRule(AmoebaVdwForce self) -> std::string const & | |
def | setEpsilonCombiningRule |
setEpsilonCombiningRule(AmoebaVdwForce self, std::string const & epsilonCombiningRule) | |
def | getEpsilonCombiningRule |
getEpsilonCombiningRule(AmoebaVdwForce self) -> std::string const & | |
def | getUseDispersionCorrection |
getUseDispersionCorrection(AmoebaVdwForce self) -> bool | |
def | setUseDispersionCorrection |
setUseDispersionCorrection(AmoebaVdwForce self, bool useCorrection) | |
def | setParticleExclusions |
setParticleExclusions(AmoebaVdwForce self, int particleIndex, vectori exclusions) | |
def | getParticleExclusions |
getParticleExclusions(AmoebaVdwForce self, int particleIndex) | |
def | setCutoff |
setCutoff(AmoebaVdwForce self, double cutoff) | |
def | getCutoff |
getCutoff(AmoebaVdwForce self) -> double | |
def | getNonbondedMethod |
getNonbondedMethod(AmoebaVdwForce self) -> OpenMM::AmoebaVdwForce::NonbondedMethod | |
def | setNonbondedMethod |
setNonbondedMethod(AmoebaVdwForce self, OpenMM::AmoebaVdwForce::NonbondedMethod method) | |
def | updateParametersInContext |
updateParametersInContext(AmoebaVdwForce self, Context context) | |
def | __init__ |
init(OpenMM::AmoebaVdwForce self) -> AmoebaVdwForce init(OpenMM::AmoebaVdwForce self, AmoebaVdwForce other) -> AmoebaVdwForce | |
def | __del__ |
del(OpenMM::AmoebaVdwForce 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.AmoebaVdwForce_NoCutoff | |
CutoffPeriodic _openmm.AmoebaVdwForce_CutoffPeriodic | |
This class implements a buffered 14-7 potential used to model van der Waals forces.
To use it, create an AmoebaVdwForce object then call addParticle() once for each particle. 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().
A unique feature of this class is that the interaction site for a particle does not need to be exactly at the particle's location. Instead, it can be placed a fraction of the distance from that particle to another one. This is typically done for hydrogens to place the interaction site slightly closer to the parent atom. The fraction is known as the "reduction factor", since it reduces the distance from the parent atom to the interaction site.
def __init__ | ( | self, | |
args | |||
) |
init(OpenMM::AmoebaVdwForce self) -> AmoebaVdwForce init(OpenMM::AmoebaVdwForce self, AmoebaVdwForce other) -> AmoebaVdwForce
Create an Amoeba VdwForce.
def __del__ | ( | self | ) |
del(OpenMM::AmoebaVdwForce self)
def addParticle | ( | self, | |
args | |||
) |
addParticle(AmoebaVdwForce self, int parentIndex, double sigma, double epsilon, double reductionFactor) -> int
Add the force field parameters for a vdw particle.
parentIndex | the index of the parent particle |
sigma | vdw sigma |
epsilon | vdw epsilon |
reductionFactor | the fraction of the distance along the line from the parent particle to this particle at which the interaction site should be placed |
def getCutoff | ( | self | ) |
getCutoff(AmoebaVdwForce self) -> double
Get the cutoff distance.
def getEpsilonCombiningRule | ( | self | ) |
getEpsilonCombiningRule(AmoebaVdwForce self) -> std::string const &
Get epsilon combining rule
def getNonbondedMethod | ( | self | ) |
getNonbondedMethod(AmoebaVdwForce self) -> OpenMM::AmoebaVdwForce::NonbondedMethod
Get the method used for handling long range nonbonded interactions.
def getNumParticles | ( | self | ) |
getNumParticles(AmoebaVdwForce self) -> int
Get the number of particles
def getParticleExclusions | ( | self, | |
args | |||
) |
getParticleExclusions(AmoebaVdwForce self, int particleIndex)
Get exclusions for specified particle
particleIndex | particle index |
exclusions | vector of exclusions |
def getParticleParameters | ( | self, | |
args | |||
) |
getParticleParameters(AmoebaVdwForce self, int particleIndex)
Get the force field parameters for a vdw particle.
particleIndex | the particle index |
parentIndex | the index of the parent particle |
sigma | vdw sigma |
epsilon | vdw epsilon |
reductionFactor | the fraction of the distance along the line from the parent particle to this particle at which the interaction site should be placed |
def getSigmaCombiningRule | ( | self | ) |
getSigmaCombiningRule(AmoebaVdwForce self) -> std::string const &
Get sigma combining rule
def getUseDispersionCorrection | ( | self | ) |
getUseDispersionCorrection(AmoebaVdwForce self) -> bool
Get whether to add a contribution to the energy that approximately represents the effect of VdW interactions beyond the cutoff distance. The energy depends on the volume of the periodic box, and is only applicable when periodic boundary conditions are used. When running simulations at constant pressure, adding this contribution can improve the quality of results.
def setCutoff | ( | self, | |
args | |||
) |
setCutoff(AmoebaVdwForce self, double cutoff)
Set the cutoff distance.
def setEpsilonCombiningRule | ( | self, | |
args | |||
) |
setEpsilonCombiningRule(AmoebaVdwForce self, std::string const & epsilonCombiningRule)
Set epsilon combining rule
epsilonCombiningRule | epsilon combining rule: 'ARITHMETIC', 'GEOMETRIC'. 'HARMONIC', 'HHG' |
def setNonbondedMethod | ( | self, | |
args | |||
) |
setNonbondedMethod(AmoebaVdwForce self, OpenMM::AmoebaVdwForce::NonbondedMethod method)
Set the method used for handling long range nonbonded interactions.
def setParticleExclusions | ( | self, | |
args | |||
) |
setParticleExclusions(AmoebaVdwForce self, int particleIndex, vectori exclusions)
Set exclusions for specified particle
particleIndex | particle index |
exclusions | vector of exclusions |
def setParticleParameters | ( | self, | |
args | |||
) |
setParticleParameters(AmoebaVdwForce self, int particleIndex, int parentIndex, double sigma, double epsilon, double reductionFactor)
Set the force field parameters for a vdw particle.
particleIndex | the particle index |
parentIndex | the index of the parent particle |
sigma | vdw sigma |
epsilon | vdw epsilon |
reductionFactor | the fraction of the distance along the line from the parent particle to this particle at which the interaction site should be placed |
def setSigmaCombiningRule | ( | self, | |
args | |||
) |
setSigmaCombiningRule(AmoebaVdwForce self, std::string const & sigmaCombiningRule)
Set sigma combining rule
sigmaCombiningRule | sigma combining rule: 'ARITHMETIC', 'GEOMETRIC'. 'CUBIC-MEAN' |
def setUseDispersionCorrection | ( | self, | |
args | |||
) |
setUseDispersionCorrection(AmoebaVdwForce self, bool useCorrection)
Set whether to add a contribution to the energy that approximately represents the effect of VdW interactions beyond the cutoff distance. The energy depends on the volume of the periodic box, and is only applicable when periodic boundary conditions are used. When running simulations at constant pressure, adding this contribution can improve the quality of results.
def updateParametersInContext | ( | self, | |
args | |||
) |
updateParametersInContext(AmoebaVdwForce self, Context context)
Update the per-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.
|
static |
|
static |
this |