AmoebaMultipoleForce Class Reference

This class implements the Amoeba multipole interaction To use it, create a MultipoleForce object then call addMultipole() once for each atom. More...

Inheritance diagram for AmoebaMultipoleForce:
Force

List of all members.

Public Member Functions

def getNumMultipoles
 getNumMultipoles(self) -> int
def getNonbondedMethod
 getNonbondedMethod(self) -> AmoebaNonbondedMethod
def setNonbondedMethod
 setNonbondedMethod(self, AmoebaNonbondedMethod method)
def getPolarizationType
 getPolarizationType(self) -> AmoebaPolarizationType
def setPolarizationType
 setPolarizationType(self, AmoebaPolarizationType type)
def getCutoffDistance
 getCutoffDistance(self) -> double
def setCutoffDistance
 setCutoffDistance(self, double distance)
def getAEwald
 getAEwald(self) -> double
def setAEwald
 setAEwald(self, double aewald)
def getPmeBSplineOrder
 getPmeBSplineOrder(self) -> int
def getPmeGridDimensions
 getPmeGridDimensions(self)
def setPmeGridDimensions
 setPmeGridDimensions(self, vectori gridDimension)
def addParticle
 addParticle(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(self, int index)
def setMultipoleParameters
 setMultipoleParameters(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(self, int index, CovalentType typeId, vectori covalentAtoms)
def getCovalentMap
 getCovalentMap(self, int index, CovalentType typeId)
def getCovalentMaps
 getCovalentMaps(self, int index)
def getMutualInducedMaxIterations
 getMutualInducedMaxIterations(self) -> int
def setMutualInducedMaxIterations
 setMutualInducedMaxIterations(self, int inputMutualInducedMaxIterations)
def getMutualInducedTargetEpsilon
 getMutualInducedTargetEpsilon(self) -> double
def setMutualInducedTargetEpsilon
 setMutualInducedTargetEpsilon(self, double inputMutualInducedTargetEpsilon)
def getEwaldErrorTolerance
 getEwaldErrorTolerance(self) -> double
def setEwaldErrorTolerance
 setEwaldErrorTolerance(self, double tol)
def __init__
 __init__(self) -> AmoebaMultipoleForce __init__(self, AmoebaMultipoleForce other) -> AmoebaMultipoleForce
def __del__
 __del__(self)

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
 SOR = _openmm.AmoebaMultipoleForce_SOR
 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 a MultipoleForce object then call addMultipole() once for each atom.

After a entry has been added, you can modify its force field parameters by calling setMultipoleParameters().


Member Function Documentation

def __del__ (   self  ) 

__del__(self)

Reimplemented from Force.

def __init__ (   self,
  args 
)

__init__(self) -> AmoebaMultipoleForce __init__(self, AmoebaMultipoleForce other) -> AmoebaMultipoleForce

Create a Amoeba MultipoleForce.

def addParticle (   self,
  args 
)

addParticle(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:
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 ( ZThenX, Bisector )
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
def getAEwald (   self  ) 

getAEwald(self) -> double

Get the aEwald parameter

def getCovalentMap (   self,
  args 
)

getCovalentMap(self, int index, CovalentType typeId)

Get the CovalentMap for an atom

Parameters:
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
def getCovalentMaps (   self,
  args 
)

getCovalentMaps(self, int index)

Get the CovalentMap for an atom

Parameters:
index the index of the atom for which to set parameters
covalentLists output vector of covalent lists of atoms
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 getEwaldErrorTolerance (   self  ) 

getEwaldErrorTolerance(self) -> double

Get the scaling distance cutoff (nm)

Parameters:
scaling distance cutoff Get the electric constant
the electric constant 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 reciprocal space cutoff and separation 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.
def getMultipoleParameters (   self,
  args 
)

getMultipoleParameters(self, int index)

Get the multipole parameters for a particle.

Parameters:
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 ( ZThenX, Bisector )
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
def getMutualInducedMaxIterations (   self  ) 

getMutualInducedMaxIterations(self) -> int

Get the iteration method to be used for calculating the mutual induced dipoles

Parameters:
iteration method to be used for calculating the mutual induced dipole 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) -> AmoebaNonbondedMethod

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 parameter

def getPmeGridDimensions (   self  ) 

getPmeGridDimensions(self)

Set the B-spline order parameter

Parameters:
the B-spline order parameter Get the PME grid dimensions
def getPolarizationType (   self  ) 

getPolarizationType(self) -> AmoebaPolarizationType

Get polarization type

def setAEwald (   self,
  args 
)

setAEwald(self, double aewald)

Set the aEwald parameter

Parameters:
Ewald parameter
def setCovalentMap (   self,
  args 
)

setCovalentMap(self, int index, CovalentType typeId, vectori covalentAtoms)

Set the CovalentMap for an atom

Parameters:
index the index of the atom for which to set parameters
typeId CovalentTypes type
covalentAtoms vector of covalent atoms associated w/ the specfied CovalentType
def setCutoffDistance (   self,
  args 
)

setCutoffDistance(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:
distance the cutoff distance, measured in nm
def setEwaldErrorTolerance (   self,
  args 
)

setEwaldErrorTolerance(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 reciprocal space cutoff and separation 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.

def setMultipoleParameters (   self,
  args 
)

setMultipoleParameters(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:
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 ( ZThenX, Bisector )
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
def setMutualInducedMaxIterations (   self,
  args 
)

setMutualInducedMaxIterations(self, int inputMutualInducedMaxIterations)

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

Parameters:
max number of iterations
def setMutualInducedTargetEpsilon (   self,
  args 
)

setMutualInducedTargetEpsilon(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:
target epsilon
def setNonbondedMethod (   self,
  args 
)

setNonbondedMethod(self, AmoebaNonbondedMethod method)

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

def setPmeGridDimensions (   self,
  args 
)

setPmeGridDimensions(self, vectori gridDimension)

Set the PME grid dimensions

Parameters:
the PME grid dimensions
def setPolarizationType (   self,
  args 
)

setPolarizationType(self, AmoebaPolarizationType type)

Set the polarization type


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]
SOR = _openmm.AmoebaMultipoleForce_SOR [static]
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:

Generated by  doxygen 1.6.2