OpenMM
|
This class implements the Amoeba multipole interaction. More...
Public Member Functions | |
def | getNumMultipoles |
getNumMultipoles(self) -> int | |
def | getNonbondedMethod |
getNonbondedMethod(self) -> OpenMM::AmoebaMultipoleForce::NonbondedMethod | |
def | setNonbondedMethod |
Set the method used for handling long-range nonbonded interactions. | |
def | getPolarizationType |
getPolarizationType(self) -> OpenMM::AmoebaMultipoleForce::PolarizationType | |
def | setPolarizationType |
Set the polarization type. | |
def | getCutoffDistance |
getCutoffDistance(self) -> double | |
def | setCutoffDistance |
Set the cutoff distance (in nm) being used for nonbonded interactions. | |
def | getAEwald |
getAEwald(self) -> double | |
def | setAEwald |
Set the Ewald alpha parameter. | |
def | getPmeBSplineOrder |
getPmeBSplineOrder(self) -> int | |
def | getPmeGridDimensions |
Get the PME grid dimensions. | |
def | setPmeGridDimensions |
Set the PME grid dimensions. | |
def | getPMEParametersInContext |
Get the parameters being used for PME in a particular Context. | |
def | addMultipole |
addMultipole(self, charge, molecularDipole, molecularQuadrupole, axisType, multipoleAtomZ, multipoleAtomX, multipoleAtomY, thole, dampingFactor, polarity) -> int | |
def | getMultipoleParameters |
Get the multipole parameters for a particle. | |
def | setMultipoleParameters |
Set the multipole parameters for a particle. | |
def | setCovalentMap |
Set the CovalentMap for an atom. | |
def | getCovalentMap |
Get the CovalentMap for an atom. | |
def | getCovalentMaps |
Get the CovalentMap for an atom. | |
def | getMutualInducedMaxIterations |
getMutualInducedMaxIterations(self) -> int | |
def | setMutualInducedMaxIterations |
Set the max number of iterations to be used in calculating the mutual induced dipoles. | |
def | getMutualInducedTargetEpsilon |
getMutualInducedTargetEpsilon(self) -> double | |
def | setMutualInducedTargetEpsilon |
Set the target epsilon to be used to test for convergence of iterative method used in calculating the mutual induced dipoles. | |
def | getEwaldErrorTolerance |
getEwaldErrorTolerance(self) -> double | |
def | setEwaldErrorTolerance |
Get the error tolerance for Ewald summation. | |
def | getInducedDipoles |
Get the induced dipole moments of all particles. | |
def | getElectrostaticPotential |
Get the electrostatic potential. | |
def | getSystemMultipoleMoments |
Get the system multipole moments. | |
def | updateParametersInContext |
Update the multipole parameters in a Context to match those stored in this Force object. | |
def | usesPeriodicBoundaryConditions |
usesPeriodicBoundaryConditions(self) -> bool | |
def | __init__ |
__init__(self) -> AmoebaMultipoleForce __init__(self, other) -> AmoebaMultipoleForce | |
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__(self) -> AmoebaMultipoleForce __init__(self, other) -> AmoebaMultipoleForce
Create an AmoebaMultipoleForce.
def addMultipole | ( | self, | |
charge, | |||
molecularDipole, | |||
molecularQuadrupole, | |||
axisType, | |||
multipoleAtomZ, | |||
multipoleAtomX, | |||
multipoleAtomY, | |||
thole, | |||
dampingFactor, | |||
polarity | |||
) |
addMultipole(self, charge, molecularDipole, molecularQuadrupole, axisType, multipoleAtomZ, multipoleAtomX, multipoleAtomY, thole, dampingFactor, polarity) -> int
Add multipole-related info for a particle
charge | (double) the particle's charge |
molecularDipole | (vector< double >) the particle's molecular dipole (vector of size 3) |
molecularQuadrupole | (vector< double >) the particle's molecular quadrupole (vector of size 9) |
axisType | (int) the particle's axis type |
multipoleAtomZ | (int) index of first atom used in constructing lab<->molecular frames |
multipoleAtomX | (int) index of second atom used in constructing lab<->molecular frames |
multipoleAtomY | (int) index of second atom used in constructing lab<->molecular frames |
thole | (double) Thole parameter |
dampingFactor | (double) dampingFactor parameter |
polarity | (double) polarity parameter |
def getAEwald | ( | self | ) |
getAEwald(self) -> double
Get the Ewald alpha parameter. If this is 0 (the default), a value is chosen automatically based on the Ewald error tolerance.
def getCovalentMap | ( | self, | |
index, | |||
typeId | |||
) |
Get the CovalentMap for an atom.
index | (int) the index of the atom for which to set parameters |
typeId | (CovalentType) CovalentTypes type |
def getCovalentMaps | ( | self, | |
index | |||
) |
Get the CovalentMap for an atom.
index | (int) the index of the atom for which to set parameters |
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 getElectrostaticPotential | ( | self, | |
inputGrid, | |||
context | |||
) |
Get the electrostatic potential.
inputGrid | (vector< Vec3 >) input grid points over which the potential is to be evaluated |
context | (Context) context |
def getEwaldErrorTolerance | ( | self | ) |
getEwaldErrorTolerance(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.
def getInducedDipoles | ( | self, | |
context | |||
) |
def getMultipoleParameters | ( | self, | |
index | |||
) |
Get the multipole parameters for a particle.
index | (int) the index of the atom for which to get parameters |
def getMutualInducedMaxIterations | ( | self | ) |
getMutualInducedMaxIterations(self) -> int
Get the max number of iterations to be used in calculating the mutual induced dipoles
def getMutualInducedTargetEpsilon | ( | self | ) |
getMutualInducedTargetEpsilon(self) -> double
Get the target epsilon to be used to test for convergence of iterative method used in calculating the mutual induced dipoles
def getNonbondedMethod | ( | self | ) |
getNonbondedMethod(self) -> OpenMM::AmoebaMultipoleForce::NonbondedMethod
Get the method used for handling long-range nonbonded interactions.
def getNumMultipoles | ( | self | ) |
getNumMultipoles(self) -> int
Get the number of particles in the potential function
def getPmeBSplineOrder | ( | self | ) |
getPmeBSplineOrder(self) -> int
Get the B-spline order to use for PME charge spreading
def getPmeGridDimensions | ( | 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.
def getPMEParametersInContext | ( | self, | |
context | |||
) |
Get the parameters being used for PME in a particular Context.
Because some platforms have restrictions on the allowed grid sizes, the values that are actually used may be slightly different from those specified with setPmeGridDimensions(), or the standard values calculated based on the Ewald error tolerance. See the manual for details.
def getPolarizationType | ( | self | ) |
getPolarizationType(self) -> OpenMM::AmoebaMultipoleForce::PolarizationType
Get polarization type
def getSystemMultipoleMoments | ( | self, | |
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) context |
def setAEwald | ( | self, | |
aewald | |||
) |
Set the Ewald alpha parameter.
If this is 0 (the default), a value is chosen automatically based on the Ewald error tolerance.
aewald | (double) alpha parameter |
def setCovalentMap | ( | self, | |
index, | |||
typeId, | |||
covalentAtoms | |||
) |
Set the CovalentMap for an atom.
index | (int) the index of the atom for which to set parameters |
typeId | (CovalentType) CovalentTypes type |
covalentAtoms | (vector< int >) vector of covalent atoms associated w/ the specfied CovalentType |
def setCutoffDistance | ( | self, | |
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 | (double) the cutoff distance, measured in nm |
def setEwaldErrorTolerance | ( | self, | |
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.
def setMultipoleParameters | ( | self, | |
index, | |||
charge, | |||
molecularDipole, | |||
molecularQuadrupole, | |||
axisType, | |||
multipoleAtomZ, | |||
multipoleAtomX, | |||
multipoleAtomY, | |||
thole, | |||
dampingFactor, | |||
polarity | |||
) |
Set the multipole parameters for a particle.
index | (int) the index of the atom for which to set parameters |
charge | (double) the particle's charge |
molecularDipole | (vector< double >) the particle's molecular dipole (vector of size 3) |
molecularQuadrupole | (vector< double >) the particle's molecular quadrupole (vector of size 9) |
axisType | (int) the particle's axis type |
multipoleAtomZ | (int) index of first atom used in constructing lab<->molecular frames |
multipoleAtomX | (int) index of second atom used in constructing lab<->molecular frames |
multipoleAtomY | (int) index of second atom used in constructing lab<->molecular frames |
thole | (double) thole parameter |
dampingFactor | (double) damping factor parameter |
polarity | (double) polarity parameter |
def setMutualInducedMaxIterations | ( | self, | |
inputMutualInducedMaxIterations | |||
) |
Set the max number of iterations to be used in calculating the mutual induced dipoles.
inputMutualInducedMaxIterations | (int) number of iterations |
def setMutualInducedTargetEpsilon | ( | self, | |
inputMutualInducedTargetEpsilon | |||
) |
Set the target epsilon to be used to test for convergence of iterative method used in calculating the mutual induced dipoles.
inputMutualInducedTargetEpsilon | (double) target epsilon |
def setNonbondedMethod | ( | self, | |
method | |||
) |
Set the method used for handling long-range nonbonded interactions.
def setPmeGridDimensions | ( | self, | |
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.
gridDimension | (vector< int >) the PME grid dimensions |
def setPolarizationType | ( | self, | |
type | |||
) |
Set the polarization type.
def updateParametersInContext | ( | self, | |
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 updateParametersInContext() 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.
def usesPeriodicBoundaryConditions | ( | self | ) |
usesPeriodicBoundaryConditions(self) -> bool
Returns whether or not this force makes use of periodic boundary conditions.
Reimplemented from Force.
Bisector = _openmm.AmoebaMultipoleForce_Bisector [static] |
Covalent12 = _openmm.AmoebaMultipoleForce_Covalent12 [static] |
Covalent13 = _openmm.AmoebaMultipoleForce_Covalent13 [static] |
Covalent14 = _openmm.AmoebaMultipoleForce_Covalent14 [static] |
Covalent15 = _openmm.AmoebaMultipoleForce_Covalent15 [static] |
CovalentEnd = _openmm.AmoebaMultipoleForce_CovalentEnd [static] |
Direct = _openmm.AmoebaMultipoleForce_Direct [static] |
LastAxisTypeIndex = _openmm.AmoebaMultipoleForce_LastAxisTypeIndex [static] |
Mutual = _openmm.AmoebaMultipoleForce_Mutual [static] |
NoAxisType = _openmm.AmoebaMultipoleForce_NoAxisType [static] |
NoCutoff = _openmm.AmoebaMultipoleForce_NoCutoff [static] |
PME = _openmm.AmoebaMultipoleForce_PME [static] |
PolarizationCovalent11 = _openmm.AmoebaMultipoleForce_PolarizationCovalent11 [static] |
PolarizationCovalent12 = _openmm.AmoebaMultipoleForce_PolarizationCovalent12 [static] |
PolarizationCovalent13 = _openmm.AmoebaMultipoleForce_PolarizationCovalent13 [static] |
PolarizationCovalent14 = _openmm.AmoebaMultipoleForce_PolarizationCovalent14 [static] |
ThreeFold = _openmm.AmoebaMultipoleForce_ThreeFold [static] |
ZBisect = _openmm.AmoebaMultipoleForce_ZBisect [static] |
ZOnly = _openmm.AmoebaMultipoleForce_ZOnly [static] |
ZThenX = _openmm.AmoebaMultipoleForce_ZThenX [static] |