OpenMM
|
This class implements the Amoeba multipole interaction. More...
Public Member Functions | |
def | getNumMultipoles |
getNumMultipoles(AmoebaMultipoleForce self) -> int More... | |
def | getNonbondedMethod |
getNonbondedMethod(AmoebaMultipoleForce self) -> OpenMM::AmoebaMultipoleForce::NonbondedMethod More... | |
def | setNonbondedMethod |
setNonbondedMethod(AmoebaMultipoleForce self, OpenMM::AmoebaMultipoleForce::NonbondedMethod method) More... | |
def | getPolarizationType |
getPolarizationType(AmoebaMultipoleForce self) -> OpenMM::AmoebaMultipoleForce::PolarizationType More... | |
def | setPolarizationType |
setPolarizationType(AmoebaMultipoleForce self, OpenMM::AmoebaMultipoleForce::PolarizationType type) More... | |
def | getCutoffDistance |
getCutoffDistance(AmoebaMultipoleForce self) -> double More... | |
def | setCutoffDistance |
setCutoffDistance(AmoebaMultipoleForce self, double distance) More... | |
def | getAEwald |
getAEwald(AmoebaMultipoleForce self) -> double More... | |
def | setAEwald |
setAEwald(AmoebaMultipoleForce self, double aewald) More... | |
def | getPmeBSplineOrder |
getPmeBSplineOrder(AmoebaMultipoleForce self) -> int More... | |
def | getPmeGridDimensions |
getPmeGridDimensions(AmoebaMultipoleForce self) More... | |
def | setPmeGridDimensions |
setPmeGridDimensions(AmoebaMultipoleForce self, vectori gridDimension) More... | |
def | addMultipole |
addMultipole(AmoebaMultipoleForce self, double charge, vectord molecularDipole, vectord molecularQuadrupole, int axisType, int multipoleAtomZ, int multipoleAtomX, int multipoleAtomY, double thole, double dampingFactor, double polarity) -> int More... | |
def | getMultipoleParameters |
getMultipoleParameters(AmoebaMultipoleForce self, int index) More... | |
def | setMultipoleParameters |
setMultipoleParameters(AmoebaMultipoleForce self, int index, double charge, vectord molecularDipole, vectord molecularQuadrupole, int axisType, int multipoleAtomZ, int multipoleAtomX, int multipoleAtomY, double thole, double dampingFactor, double polarity) More... | |
def | setCovalentMap |
setCovalentMap(AmoebaMultipoleForce self, int index, OpenMM::AmoebaMultipoleForce::CovalentType typeId, vectori covalentAtoms) More... | |
def | getCovalentMap |
getCovalentMap(AmoebaMultipoleForce self, int index, OpenMM::AmoebaMultipoleForce::CovalentType typeId) More... | |
def | getCovalentMaps |
getCovalentMaps(AmoebaMultipoleForce self, int index) More... | |
def | getMutualInducedMaxIterations |
getMutualInducedMaxIterations(AmoebaMultipoleForce self) -> int More... | |
def | setMutualInducedMaxIterations |
setMutualInducedMaxIterations(AmoebaMultipoleForce self, int inputMutualInducedMaxIterations) More... | |
def | getMutualInducedTargetEpsilon |
getMutualInducedTargetEpsilon(AmoebaMultipoleForce self) -> double More... | |
def | setMutualInducedTargetEpsilon |
setMutualInducedTargetEpsilon(AmoebaMultipoleForce self, double inputMutualInducedTargetEpsilon) More... | |
def | getEwaldErrorTolerance |
getEwaldErrorTolerance(AmoebaMultipoleForce self) -> double More... | |
def | setEwaldErrorTolerance |
setEwaldErrorTolerance(AmoebaMultipoleForce self, double tol) More... | |
def | getInducedDipoles |
getInducedDipoles(AmoebaMultipoleForce self, Context context) More... | |
def | getElectrostaticPotential |
getElectrostaticPotential(AmoebaMultipoleForce self, std::vector< Vec3,std::allocator< Vec3 > > const & inputGrid, Context context) More... | |
def | getSystemMultipoleMoments |
getSystemMultipoleMoments(AmoebaMultipoleForce self, Context context) More... | |
def | updateParametersInContext |
updateParametersInContext(AmoebaMultipoleForce self, Context context) More... | |
def | __init__ |
init(OpenMM::AmoebaMultipoleForce self) -> AmoebaMultipoleForce init(OpenMM::AmoebaMultipoleForce self, AmoebaMultipoleForce other) -> AmoebaMultipoleForce More... | |
def | __del__ |
del(OpenMM::AmoebaMultipoleForce self) More... | |
Public Member Functions inherited from Force | |
def | __init__ |
def | __del__ |
del(OpenMM::Force self) More... | |
def | getForceGroup |
getForceGroup(Force self) -> int More... | |
def | setForceGroup |
setForceGroup(Force self, int group) More... | |
def | __copy__ |
def | __deepcopy__ |
Public Attributes | |
this | |
Static Public Attributes | |
NoCutoff = _openmm.AmoebaMultipoleForce_NoCutoff | |
PME = _openmm.AmoebaMultipoleForce_PME | |
Mutual = _openmm.AmoebaMultipoleForce_Mutual | |
Direct = _openmm.AmoebaMultipoleForce_Direct | |
ZThenX = _openmm.AmoebaMultipoleForce_ZThenX | |
Bisector = _openmm.AmoebaMultipoleForce_Bisector | |
ZBisect = _openmm.AmoebaMultipoleForce_ZBisect | |
ThreeFold = _openmm.AmoebaMultipoleForce_ThreeFold | |
ZOnly = _openmm.AmoebaMultipoleForce_ZOnly | |
NoAxisType = _openmm.AmoebaMultipoleForce_NoAxisType | |
LastAxisTypeIndex = _openmm.AmoebaMultipoleForce_LastAxisTypeIndex | |
Covalent12 = _openmm.AmoebaMultipoleForce_Covalent12 | |
Covalent13 = _openmm.AmoebaMultipoleForce_Covalent13 | |
Covalent14 = _openmm.AmoebaMultipoleForce_Covalent14 | |
Covalent15 = _openmm.AmoebaMultipoleForce_Covalent15 | |
PolarizationCovalent11 = _openmm.AmoebaMultipoleForce_PolarizationCovalent11 | |
PolarizationCovalent12 = _openmm.AmoebaMultipoleForce_PolarizationCovalent12 | |
PolarizationCovalent13 = _openmm.AmoebaMultipoleForce_PolarizationCovalent13 | |
PolarizationCovalent14 = _openmm.AmoebaMultipoleForce_PolarizationCovalent14 | |
CovalentEnd = _openmm.AmoebaMultipoleForce_CovalentEnd | |
This class implements the Amoeba multipole interaction.
To use it, create an AmoebaMultipoleForce object then call addMultipole() once for each atom. After an entry has been added, you can modify its force field parameters by calling setMultipoleParameters(). This will have no effect on Contexts that already exist unless you call updateParametersInContext().
def __init__ | ( | self, | |
args | |||
) |
init(OpenMM::AmoebaMultipoleForce self) -> AmoebaMultipoleForce init(OpenMM::AmoebaMultipoleForce self, AmoebaMultipoleForce other) -> AmoebaMultipoleForce
Create an AmoebaMultipoleForce.
References simtk.openmm.openmm.stripUnits().
def __del__ | ( | self | ) |
del(OpenMM::AmoebaMultipoleForce self)
References simtk.openmm.openmm.stripUnits().
def addMultipole | ( | self, | |
args | |||
) |
addMultipole(AmoebaMultipoleForce self, double charge, vectord molecularDipole, vectord molecularQuadrupole, int axisType, int multipoleAtomZ, int multipoleAtomX, int multipoleAtomY, double thole, double dampingFactor, double polarity) -> int
Add multipole-related info for a particle
charge | the particle's charge |
molecularDipole | the particle's molecular dipole (vector of size 3) |
molecularQuadrupole | the particle's molecular quadrupole (vector of size 9) |
axisType | the particle's axis type |
multipoleAtomZ | index of first atom used in constructing lab<->molecular frames |
multipoleAtomX | index of second atom used in constructing lab<->molecular frames |
multipoleAtomY | index of second atom used in constructing lab<->molecular frames |
thole | Thole parameter |
dampingFactor | dampingFactor parameter |
polarity | polarity parameter |
References simtk.openmm.openmm.stripUnits().
def getAEwald | ( | self | ) |
getAEwald(AmoebaMultipoleForce self) -> double
Get the Ewald alpha parameter. If this is 0 (the default), a value is chosen automatically based on the Ewald error tolerance.
References simtk.openmm.openmm.stripUnits().
def getCovalentMap | ( | self, | |
args | |||
) |
getCovalentMap(AmoebaMultipoleForce self, int index, OpenMM::AmoebaMultipoleForce::CovalentType typeId)
Get the CovalentMap for an atom
index | the index of the atom for which to set parameters |
typeId | CovalentTypes type |
covalentAtoms | output vector of covalent atoms associated w/ the specfied CovalentType |
References simtk.openmm.openmm.stripUnits().
def getCovalentMaps | ( | self, | |
args | |||
) |
getCovalentMaps(AmoebaMultipoleForce self, int index)
Get the CovalentMap for an atom
index | the index of the atom for which to set parameters |
covalentLists | output vector of covalent lists of atoms |
References simtk.openmm.openmm.stripUnits().
def getCutoffDistance | ( | self | ) |
getCutoffDistance(AmoebaMultipoleForce 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.
References simtk.openmm.openmm.stripUnits().
def getElectrostaticPotential | ( | self, | |
args | |||
) |
getElectrostaticPotential(AmoebaMultipoleForce self, std::vector< Vec3,std::allocator< Vec3 > > const & inputGrid, Context context)
Get the electrostatic potential.
inputGrid | input grid points over which the potential is to be evaluated |
context | context |
outputElectrostaticPotential | output potential |
References simtk.openmm.openmm.stripUnits().
def getEwaldErrorTolerance | ( | self | ) |
getEwaldErrorTolerance(AmoebaMultipoleForce self) -> double
Get the error tolerance for Ewald summation. This corresponds to the fractional error in the forces which is acceptable. This value is used to select the grid dimensions and separation (alpha) parameter so that the average error level will be less than the tolerance. There is not a rigorous guarantee that all forces on all atoms will be less than the tolerance, however.
This can be overridden by explicitly setting an alpha parameter and grid dimensions to use.
References simtk.openmm.openmm.stripUnits().
def getInducedDipoles | ( | self, | |
args | |||
) |
getInducedDipoles(AmoebaMultipoleForce self, Context context)
Get the induced dipole moments of all particles.
context | the Context for which to get the induced dipoles |
dipoles | the induced dipole moment of particle i is stored into the i'th element |
References simtk.openmm.openmm.stripUnits().
def getMultipoleParameters | ( | self, | |
args | |||
) |
getMultipoleParameters(AmoebaMultipoleForce self, int index)
Get the multipole parameters for a particle.
index | the index of the atom for which to get parameters |
charge | the particle's charge |
molecularDipole | the particle's molecular dipole (vector of size 3) |
molecularQuadrupole | the particle's molecular quadrupole (vector of size 9) |
axisType | the particle's axis type |
multipoleAtomZ | index of first atom used in constructing lab<->molecular frames |
multipoleAtomX | index of second atom used in constructing lab<->molecular frames |
multipoleAtomY | index of second atom used in constructing lab<->molecular frames |
thole | Thole parameter |
dampingFactor | dampingFactor parameter |
polarity | polarity parameter |
References simtk.openmm.openmm.stripUnits().
def getMutualInducedMaxIterations | ( | self | ) |
getMutualInducedMaxIterations(AmoebaMultipoleForce self) -> int
Get the max number of iterations to be used in calculating the mutual induced dipoles
References simtk.openmm.openmm.stripUnits().
def getMutualInducedTargetEpsilon | ( | self | ) |
getMutualInducedTargetEpsilon(AmoebaMultipoleForce self) -> double
Get the target epsilon to be used to test for convergence of iterative method used in calculating the mutual induced dipoles
References simtk.openmm.openmm.stripUnits().
def getNonbondedMethod | ( | self | ) |
getNonbondedMethod(AmoebaMultipoleForce self) -> OpenMM::AmoebaMultipoleForce::NonbondedMethod
Get the method used for handling long-range nonbonded interactions.
References simtk.openmm.openmm.stripUnits().
def getNumMultipoles | ( | self | ) |
getNumMultipoles(AmoebaMultipoleForce self) -> int
Get the number of particles in the potential function
References simtk.openmm.openmm.stripUnits().
def getPmeBSplineOrder | ( | self | ) |
getPmeBSplineOrder(AmoebaMultipoleForce self) -> int
Get the B-spline order to use for PME charge spreading
References simtk.openmm.openmm.stripUnits().
def getPmeGridDimensions | ( | self | ) |
getPmeGridDimensions(AmoebaMultipoleForce self)
Get the PME grid dimensions. If Ewald alpha is 0 (the default), this is ignored and grid dimensions are chosen automatically based on the Ewald error tolerance.
References simtk.openmm.openmm.stripUnits().
def getPolarizationType | ( | self | ) |
getPolarizationType(AmoebaMultipoleForce self) -> OpenMM::AmoebaMultipoleForce::PolarizationType
Get polarization type
References simtk.openmm.openmm.stripUnits().
def getSystemMultipoleMoments | ( | self, | |
args | |||
) |
getSystemMultipoleMoments(AmoebaMultipoleForce self, Context context)
Get the system multipole moments.
This method is most useful for non-periodic systems. When called for a periodic system, only the lowest nonvanishing moment has a well defined value. This means that if the system has a net nonzero charge, the dipole and quadrupole moments are not well defined and should be ignored. If the net charge is zero, the dipole moment is well defined (and really represents a dipole density), but the quadrupole moment is still undefined and should be ignored.
context | context |
outputMultipoleMoments | (charge, dipole_x, dipole_y, dipole_z, quadrupole_xx, quadrupole_xy, quadrupole_xz, quadrupole_yx, quadrupole_yy, quadrupole_yz, quadrupole_zx, quadrupole_zy, quadrupole_zz) |
References simtk.openmm.openmm.stripUnits().
def setAEwald | ( | self, | |
args | |||
) |
setAEwald(AmoebaMultipoleForce self, double aewald)
Set the Ewald alpha parameter. If this is 0 (the default), a value is chosen automatically based on the Ewald error tolerance.
Ewald | alpha parameter |
References simtk.openmm.openmm.stripUnits().
def setCovalentMap | ( | self, | |
args | |||
) |
setCovalentMap(AmoebaMultipoleForce self, int index, OpenMM::AmoebaMultipoleForce::CovalentType typeId, vectori covalentAtoms)
Set the CovalentMap for an atom
index | the index of the atom for which to set parameters |
typeId | CovalentTypes type |
covalentAtoms | vector of covalent atoms associated w/ the specfied CovalentType |
References simtk.openmm.openmm.stripUnits().
def setCutoffDistance | ( | self, | |
args | |||
) |
setCutoffDistance(AmoebaMultipoleForce 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 |
References simtk.openmm.openmm.stripUnits().
def setEwaldErrorTolerance | ( | self, | |
args | |||
) |
setEwaldErrorTolerance(AmoebaMultipoleForce self, double tol)
Get the error tolerance for Ewald summation. This corresponds to the fractional error in the forces which is acceptable. This value is used to select the grid dimensions and separation (alpha) parameter so that the average error level will be less than the tolerance. There is not a rigorous guarantee that all forces on all atoms will be less than the tolerance, however.
This can be overridden by explicitly setting an alpha parameter and grid dimensions to use.
References simtk.openmm.openmm.stripUnits().
def setMultipoleParameters | ( | self, | |
args | |||
) |
setMultipoleParameters(AmoebaMultipoleForce self, int index, double charge, vectord molecularDipole, vectord molecularQuadrupole, int axisType, int multipoleAtomZ, int multipoleAtomX, int multipoleAtomY, double thole, double dampingFactor, double polarity)
Set the multipole parameters for a particle.
index | the index of the atom for which to set parameters |
charge | the particle's charge |
molecularDipole | the particle's molecular dipole (vector of size 3) |
molecularQuadrupole | the particle's molecular quadrupole (vector of size 9) |
axisType | the particle's axis type |
multipoleAtomZ | index of first atom used in constructing lab<->molecular frames |
multipoleAtomX | index of second atom used in constructing lab<->molecular frames |
multipoleAtomY | index of second atom used in constructing lab<->molecular frames |
polarity | polarity parameter |
References simtk.openmm.openmm.stripUnits().
def setMutualInducedMaxIterations | ( | self, | |
args | |||
) |
setMutualInducedMaxIterations(AmoebaMultipoleForce self, int inputMutualInducedMaxIterations)
Set the max number of iterations to be used in calculating the mutual induced dipoles
max | number of iterations |
References simtk.openmm.openmm.stripUnits().
def setMutualInducedTargetEpsilon | ( | self, | |
args | |||
) |
setMutualInducedTargetEpsilon(AmoebaMultipoleForce self, double inputMutualInducedTargetEpsilon)
Set the target epsilon to be used to test for convergence of iterative method used in calculating the mutual induced dipoles
target | epsilon |
References simtk.openmm.openmm.stripUnits().
def setNonbondedMethod | ( | self, | |
args | |||
) |
setNonbondedMethod(AmoebaMultipoleForce self, OpenMM::AmoebaMultipoleForce::NonbondedMethod method)
Set the method used for handling long-range nonbonded interactions.
References simtk.openmm.openmm.stripUnits().
def setPmeGridDimensions | ( | self, | |
args | |||
) |
setPmeGridDimensions(AmoebaMultipoleForce self, vectori gridDimension)
Set the PME grid dimensions. If Ewald alpha is 0 (the default), this is ignored and grid dimensions are chosen automatically based on the Ewald error tolerance.
the | PME grid dimensions |
References simtk.openmm.openmm.stripUnits().
def setPolarizationType | ( | self, | |
args | |||
) |
setPolarizationType(AmoebaMultipoleForce self, OpenMM::AmoebaMultipoleForce::PolarizationType type)
Set the polarization type
References simtk.openmm.openmm.stripUnits().
def updateParametersInContext | ( | self, | |
args | |||
) |
updateParametersInContext(AmoebaMultipoleForce self, Context context)
Update the multipole 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 setMultipoleParameters() to modify this object's parameters, then call updateParametersInState() to copy them over to the Context.
This method has several limitations. The only information it updates is the parameters of multipoles. 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 multipoles, only to change the parameters of existing ones.
References simtk.openmm.openmm.stripUnits().
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
this |
|
static |
|
static |
|
static |
|
static |