OpenMM
 All Classes Namespaces Functions Variables Pages
AmoebaMultipoleForce Class Reference

This class implements the Amoeba multipole interaction. More...

+ Inheritance diagram for AmoebaMultipoleForce:

Public Member Functions

def getNumMultipoles
 getNumMultipoles(AmoebaMultipoleForce self) -> int
 
def getNonbondedMethod
 getNonbondedMethod(AmoebaMultipoleForce self) -> OpenMM::AmoebaMultipoleForce::NonbondedMethod
 
def setNonbondedMethod
 setNonbondedMethod(AmoebaMultipoleForce self, OpenMM::AmoebaMultipoleForce::NonbondedMethod method)
 
def getPolarizationType
 getPolarizationType(AmoebaMultipoleForce self) -> OpenMM::AmoebaMultipoleForce::PolarizationType
 
def setPolarizationType
 setPolarizationType(AmoebaMultipoleForce self, OpenMM::AmoebaMultipoleForce::PolarizationType type)
 
def getCutoffDistance
 getCutoffDistance(AmoebaMultipoleForce self) -> double
 
def setCutoffDistance
 setCutoffDistance(AmoebaMultipoleForce self, double distance)
 
def getAEwald
 getAEwald(AmoebaMultipoleForce self) -> double
 
def setAEwald
 setAEwald(AmoebaMultipoleForce self, double aewald)
 
def getPmeBSplineOrder
 getPmeBSplineOrder(AmoebaMultipoleForce self) -> int
 
def getPmeGridDimensions
 getPmeGridDimensions(AmoebaMultipoleForce self)
 
def setPmeGridDimensions
 setPmeGridDimensions(AmoebaMultipoleForce self, vectori gridDimension)
 
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
 
def getMultipoleParameters
 getMultipoleParameters(AmoebaMultipoleForce self, int index)
 
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)
 
def setCovalentMap
 setCovalentMap(AmoebaMultipoleForce self, int index, OpenMM::AmoebaMultipoleForce::CovalentType typeId, vectori covalentAtoms)
 
def getCovalentMap
 getCovalentMap(AmoebaMultipoleForce self, int index, OpenMM::AmoebaMultipoleForce::CovalentType typeId)
 
def getCovalentMaps
 getCovalentMaps(AmoebaMultipoleForce self, int index)
 
def getMutualInducedMaxIterations
 getMutualInducedMaxIterations(AmoebaMultipoleForce self) -> int
 
def setMutualInducedMaxIterations
 setMutualInducedMaxIterations(AmoebaMultipoleForce self, int inputMutualInducedMaxIterations)
 
def getMutualInducedTargetEpsilon
 getMutualInducedTargetEpsilon(AmoebaMultipoleForce self) -> double
 
def setMutualInducedTargetEpsilon
 setMutualInducedTargetEpsilon(AmoebaMultipoleForce self, double inputMutualInducedTargetEpsilon)
 
def getEwaldErrorTolerance
 getEwaldErrorTolerance(AmoebaMultipoleForce self) -> double
 
def setEwaldErrorTolerance
 setEwaldErrorTolerance(AmoebaMultipoleForce self, double tol)
 
def getElectrostaticPotential
 getElectrostaticPotential(AmoebaMultipoleForce self, std::vector< Vec3,std::allocator< Vec3 > > const & inputGrid, Context context)
 
def getSystemMultipoleMoments
 getSystemMultipoleMoments(AmoebaMultipoleForce self, Context context)
 
def updateParametersInContext
 updateParametersInContext(AmoebaMultipoleForce self, Context context)
 
def __init__
 init(OpenMM::AmoebaMultipoleForce self) -> AmoebaMultipoleForce init(OpenMM::AmoebaMultipoleForce self, AmoebaMultipoleForce other) -> AmoebaMultipoleForce
 
def __del__
 del(OpenMM::AmoebaMultipoleForce 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.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
 

Detailed Description

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().

Constructor & Destructor Documentation

def __init__ (   self,
  args 
)

init(OpenMM::AmoebaMultipoleForce self) -> AmoebaMultipoleForce init(OpenMM::AmoebaMultipoleForce self, AmoebaMultipoleForce other) -> AmoebaMultipoleForce

Create an AmoebaMultipoleForce.

def __del__ (   self)

del(OpenMM::AmoebaMultipoleForce self)

Member Function Documentation

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

Parameters
chargethe particle's charge
molecularDipolethe particle's molecular dipole (vector of size 3)
molecularQuadrupolethe particle's molecular quadrupole (vector of size 9)
axisTypethe particle's axis type
multipoleAtomZindex of first atom used in constructing lab<->molecular frames
multipoleAtomXindex of second atom used in constructing lab<->molecular frames
multipoleAtomYindex of second atom used in constructing lab<->molecular frames
tholeThole parameter
dampingFactordampingFactor parameter
polaritypolarity parameter
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.

def getCovalentMap (   self,
  args 
)

getCovalentMap(AmoebaMultipoleForce self, int index, OpenMM::AmoebaMultipoleForce::CovalentType typeId)

Get the CovalentMap for an atom

Parameters
indexthe index of the atom for which to set parameters
typeIdCovalentTypes type
covalentAtomsoutput vector of covalent atoms associated w/ the specfied CovalentType
def getCovalentMaps (   self,
  args 
)

getCovalentMaps(AmoebaMultipoleForce self, int index)

Get the CovalentMap for an atom

Parameters
indexthe index of the atom for which to set parameters
covalentListsoutput vector of covalent lists of atoms
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.

def getElectrostaticPotential (   self,
  args 
)

getElectrostaticPotential(AmoebaMultipoleForce self, std::vector< Vec3,std::allocator< Vec3 > > const & inputGrid, Context context)

Get the electrostatic potential.

Parameters
inputGridinput grid points over which the potential is to be evaluated
contextcontext
outputElectrostaticPotentialoutput potential
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.

def getMultipoleParameters (   self,
  args 
)

getMultipoleParameters(AmoebaMultipoleForce self, int index)

Get the multipole parameters for a particle.

Parameters
indexthe index of the atom for which to get parameters
chargethe particle's charge
molecularDipolethe particle's molecular dipole (vector of size 3)
molecularQuadrupolethe particle's molecular quadrupole (vector of size 9)
axisTypethe particle's axis type
multipoleAtomZindex of first atom used in constructing lab<->molecular frames
multipoleAtomXindex of second atom used in constructing lab<->molecular frames
multipoleAtomYindex of second atom used in constructing lab<->molecular frames
tholeThole parameter
dampingFactordampingFactor parameter
polaritypolarity parameter
def getMutualInducedMaxIterations (   self)

getMutualInducedMaxIterations(AmoebaMultipoleForce self) -> int

Get the max number of iterations to be used in calculating the mutual induced dipoles

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

def getNonbondedMethod (   self)

getNonbondedMethod(AmoebaMultipoleForce self) -> OpenMM::AmoebaMultipoleForce::NonbondedMethod

Get the method used for handling long-range nonbonded interactions.

def getNumMultipoles (   self)

getNumMultipoles(AmoebaMultipoleForce self) -> int

Get the number of particles in the potential function

def getPmeBSplineOrder (   self)

getPmeBSplineOrder(AmoebaMultipoleForce self) -> int

Get the B-spline order to use for PME charge spreading

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.

def getPolarizationType (   self)

getPolarizationType(AmoebaMultipoleForce self) -> OpenMM::AmoebaMultipoleForce::PolarizationType

Get polarization type

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.

Parameters
contextcontext
outputMultipoleMonents(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)
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.

Parameters
Ewaldalpha parameter
def setCovalentMap (   self,
  args 
)

setCovalentMap(AmoebaMultipoleForce self, int index, OpenMM::AmoebaMultipoleForce::CovalentType typeId, vectori covalentAtoms)

Set the CovalentMap for an atom

Parameters
indexthe index of the atom for which to set parameters
typeIdCovalentTypes type
covalentAtomsvector of covalent atoms associated w/ the specfied CovalentType
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.

Parameters
distancethe cutoff distance, measured in nm
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.

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.

Parameters
indexthe index of the atom for which to set parameters
chargethe particle's charge
molecularDipolethe particle's molecular dipole (vector of size 3)
molecularQuadrupolethe particle's molecular quadrupole (vector of size 9)
axisTypethe particle's axis type
multipoleAtomZindex of first atom used in constructing lab<->molecular frames
multipoleAtomXindex of second atom used in constructing lab<->molecular frames
multipoleAtomYindex of second atom used in constructing lab<->molecular frames
polaritypolarity parameter
def setMutualInducedMaxIterations (   self,
  args 
)

setMutualInducedMaxIterations(AmoebaMultipoleForce self, int inputMutualInducedMaxIterations)

Set the max number of iterations to be used in calculating the mutual induced dipoles

Parameters
maxnumber of iterations
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

Parameters
targetepsilon
def setNonbondedMethod (   self,
  args 
)

setNonbondedMethod(AmoebaMultipoleForce self, OpenMM::AmoebaMultipoleForce::NonbondedMethod method)

Set the method used for handling long-range nonbonded interactions.

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.

Parameters
thePME grid dimensions
def setPolarizationType (   self,
  args 
)

setPolarizationType(AmoebaMultipoleForce self, OpenMM::AmoebaMultipoleForce::PolarizationType type)

Set the polarization type

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.

Member Data Documentation

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
this
ThreeFold _openmm.AmoebaMultipoleForce_ThreeFold
static
ZBisect _openmm.AmoebaMultipoleForce_ZBisect
static
ZOnly _openmm.AmoebaMultipoleForce_ZOnly
static
ZThenX _openmm.AmoebaMultipoleForce_ZThenX
static

The documentation for this class was generated from the following file: