AmoebaMultipoleForce

class openmm.openmm.AmoebaMultipoleForce(*args)

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

__init__(self) AmoebaMultipoleForce
__init__(self, other) AmoebaMultipoleForce

Create an AmoebaMultipoleForce.

Methods

__init__(-> AmoebaMultipoleForce)

Create an AmoebaMultipoleForce.

addMultipole(self, charge, molecularDipole, ...)

Add multipole-related info for a particle

getAEwald(self)

Get the Ewald alpha parameter.

getCovalentMap(self, index, typeId)

Get the CovalentMap for an atom

getCovalentMaps(self, index)

Get the CovalentMap for an atom

getCutoffDistance(self)

Get the cutoff distance (in nm) being used for nonbonded interactions.

getElectrostaticPotential(self, inputGrid, ...)

Get the electrostatic potential.

getEwaldErrorTolerance(self)

Get the error tolerance for Ewald summation.

getExtrapolationCoefficients(self)

Get the coefficients for the mu_0, mu_1, mu_2, ..., mu_n terms in the extrapolation algorithm for induced dipoles.

getForceGroup(self)

Get the force group this Force belongs to.

getInducedDipoles(self, context)

Get the induced dipole moments of all particles.

getLabFramePermanentDipoles(self, context)

Get the fixed dipole moments of all particles in the global reference frame.

getMultipoleParameters(self, index)

Get the multipole parameters for a particle.

getMutualInducedMaxIterations(self)

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

getMutualInducedTargetEpsilon(self)

Get the target epsilon to be used to test for convergence of iterative method used in calculating the mutual induced dipoles

getName(self)

Get the name of this Force.

getNonbondedMethod(self)

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

getNumMultipoles(self)

Get the number of particles in the potential function

getPMEParameters(self)

Get the parameters to use for PME calculations.

getPMEParametersInContext(self, context)

Get the parameters being used for PME in a particular Context.

getPmeBSplineOrder(self)

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

getPmeGridDimensions(self)

Get the PME grid dimensions.

getPolarizationType(self)

Get polarization type

getSystemMultipoleMoments(self, context)

Get the system multipole moments.

getTotalDipoles(self, context)

Get the total dipole moments (fixed plus induced) of all particles.

setAEwald(self, aewald)

Set the Ewald alpha parameter.

setCovalentMap(self, index, typeId, ...)

Set the CovalentMap for an atom

setCutoffDistance(self, distance)

Set the cutoff distance (in nm) being used for nonbonded interactions.

setEwaldErrorTolerance(self, tol)

Get the error tolerance for Ewald summation.

setExtrapolationCoefficients(self, coefficients)

Set the coefficients for the mu_0, mu_1, mu_2, ..., mu_n terms in the extrapolation algorithm for induced dipoles.

setForceGroup(self, group)

Set the force group this Force belongs to.

setMultipoleParameters(self, index, charge, ...)

Set the multipole parameters for a particle.

setMutualInducedMaxIterations(self, ...)

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

setMutualInducedTargetEpsilon(self, ...)

Set the target epsilon to be used to test for convergence of iterative method used in calculating the mutual induced dipoles

setName(self, name)

Set the name of this Force.

setNonbondedMethod(self, method)

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

setPMEParameters(self, alpha, nx, ny, nz)

Set the parameters to use for PME calculations.

setPmeGridDimensions(self, gridDimension)

Set the PME grid dimensions.

setPolarizationType(self, type)

Set the polarization type

updateParametersInContext(self, context)

Update the multipole parameters in a Context to match those stored in this Force object.

usesPeriodicBoundaryConditions(self)

Returns whether or not this force makes use of periodic boundary conditions.

Attributes

Bisector

Covalent12

Covalent13

Covalent14

Covalent15

CovalentEnd

Direct

Extrapolated

LastAxisTypeIndex

Mutual

NoAxisType

NoCutoff

PME

PolarizationCovalent11

PolarizationCovalent12

PolarizationCovalent13

PolarizationCovalent14

ThreeFold

ZBisect

ZOnly

ZThenX

thisown

The membership flag

property thisown

The membership flag

getNumMultipoles(self) int

Get the number of particles in the potential function

getNonbondedMethod(self) OpenMM::AmoebaMultipoleForce::NonbondedMethod

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

setNonbondedMethod(self, method)

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

getPolarizationType(self) OpenMM::AmoebaMultipoleForce::PolarizationType

Get polarization type

setPolarizationType(self, type)

Set the polarization type

getCutoffDistance(self) float

Get the cutoff distance (in nm) being used for nonbonded interactions. If the NonbondedMethod in use is NoCutoff, this value will have no effect.

Returns:

the cutoff distance, measured in nm

Return type:

float

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.

Parameters:

distance (float) – the cutoff distance, measured in nm

getPMEParameters(self)

Get the parameters to use for PME calculations. If alpha is 0 (the default), these parameters are ignored and instead their values are chosen based on the Ewald error tolerance.

Returns:

  • alpha (float) – the separation parameter

  • nx (int) – the number of grid points along the X axis

  • ny (int) – the number of grid points along the Y axis

  • nz (int) – the number of grid points along the Z axis

setPMEParameters(self, alpha, nx, ny, nz)

Set the parameters to use for PME calculations. If alpha is 0 (the default), these parameters are ignored and instead their values are chosen based on the Ewald error tolerance.

Parameters:
  • alpha (float) – the separation parameter

  • nx (int) – the number of grid points along the X axis

  • ny (int) – the number of grid points along the Y axis

  • nz (int) – the number of grid points along the Z axis

getAEwald(self) float

Get the Ewald alpha parameter. If this is 0 (the default), a value is chosen automatically based on the Ewald error tolerance.

Deprecated

This method exists only for backward compatibility. Use getPMEParameters() instead.

Returns:

the Ewald alpha parameter

Return type:

float

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.

Deprecated

This method exists only for backward compatibility. Use setPMEParameters() instead.

Parameters:

aewald (float) – alpha parameter

getPmeBSplineOrder(self) int

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

Returns:

the B-spline order

Return type:

int

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.

Deprecated

This method exists only for backward compatibility. Use getPMEParameters() instead.

Returns:

the PME grid dimensions

Return type:

void

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.

Deprecated

This method exists only for backward compatibility. Use setPMEParameters() instead.

Parameters:

gridDimension (Sequence[int]) – the PME grid dimensions

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.

Parameters:

context (Context) – the Context for which to get the parameters

Returns:

  • alpha (float) – the separation parameter

  • nx (int) – the number of grid points along the X axis

  • ny (int) – the number of grid points along the Y axis

  • nz (int) – the number of grid points along the Z axis

addMultipole(self, charge, molecularDipole, molecularQuadrupole, axisType, multipoleAtomZ, multipoleAtomX, multipoleAtomY, thole, dampingFactor, polarity) int

Add multipole-related info for a particle

Parameters:
  • charge (float) – the particle’s charge

  • molecularDipole (Sequence[float]) – the particle’s molecular dipole (vector of size 3)

  • molecularQuadrupole (Sequence[float]) – 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 (float) – Thole parameter

  • dampingFactor (float) – dampingFactor parameter

  • polarity (float) – polarity parameter

Returns:

the index of the particle that was added

Return type:

int

getMultipoleParameters(self, index)

Get the multipole parameters for a particle.

Parameters:

index (int) – the index of the atom for which to get parameters

Returns:

  • charge (float) – the particle’s charge

  • molecularDipole (Sequence[float]) – the particle’s molecular dipole (vector of size 3)

  • molecularQuadrupole (Sequence[float]) – 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 (float) – Thole parameter

  • dampingFactor (float) – dampingFactor parameter

  • polarity (float) – polarity parameter

setMultipoleParameters(self, index, charge, molecularDipole, molecularQuadrupole, axisType, multipoleAtomZ, multipoleAtomX, multipoleAtomY, thole, dampingFactor, polarity)

Set the multipole parameters for a particle.

Parameters:
  • index (int) – the index of the atom for which to set parameters

  • charge (float) – the particle’s charge

  • molecularDipole (Sequence[float]) – the particle’s molecular dipole (vector of size 3)

  • molecularQuadrupole (Sequence[float]) – 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 (float) – thole parameter

  • dampingFactor (float) – damping factor parameter

  • polarity (float) – polarity parameter

setCovalentMap(self, index, typeId, covalentAtoms)

Set the CovalentMap for an atom

Parameters:
  • index (int) – the index of the atom for which to set parameters

  • typeId (CovalentType) – CovalentTypes type

  • covalentAtoms (Sequence[int]) – vector of covalent atoms associated w/ the specfied CovalentType

getCovalentMap(self, index, typeId)

Get the CovalentMap for an atom

Parameters:
  • index (int) – the index of the atom for which to set parameters

  • typeId (CovalentType) – CovalentTypes type

Returns:

covalentAtoms – output vector of covalent atoms associated w/ the specfied CovalentType

Return type:

Sequence[int]

getCovalentMaps(self, index)

Get the CovalentMap for an atom

Parameters:

index (int) – the index of the atom for which to set parameters

Returns:

covalentLists – output vector of covalent lists of atoms

Return type:

Sequence[Sequence[int]]

getMutualInducedMaxIterations(self) int

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

Returns:

max number of iterations

Return type:

int

setMutualInducedMaxIterations(self, inputMutualInducedMaxIterations)

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

Parameters:

inputMutualInducedMaxIterations (int) – number of iterations

getMutualInducedTargetEpsilon(self) float

Get the target epsilon to be used to test for convergence of iterative method used in calculating the mutual induced dipoles

Returns:

target epsilon

Return type:

float

setMutualInducedTargetEpsilon(self, inputMutualInducedTargetEpsilon)

Set the target epsilon to be used to test for convergence of iterative method used in calculating the mutual induced dipoles

Parameters:

inputMutualInducedTargetEpsilon (float) – target epsilon

setExtrapolationCoefficients(self, coefficients)

Set the coefficients for the mu_0, mu_1, mu_2, …, mu_n terms in the extrapolation algorithm for induced dipoles.

Parameters:

coefficients (Sequence[float]) – a vector whose mth entry specifies the coefficient for mu_m. The length of this vector determines how many iterations are performed.

getExtrapolationCoefficients(self) tuple[float, ...]

Get the coefficients for the mu_0, mu_1, mu_2, …, mu_n terms in the extrapolation algorithm for induced dipoles. In this release, the default values for the coefficients are [-0.154, 0.017, 0.658, 0.474], but be aware that those may change in a future release.

getEwaldErrorTolerance(self) float

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.

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.

getLabFramePermanentDipoles(self, context)

Get the fixed dipole moments of all particles in the global reference frame.

Parameters:

context (Context) – the Context for which to get the fixed dipoles

Returns:

dipoles – the fixed dipole moment of particle i is stored into the i’th element

Return type:

Sequence[Vec3]

getInducedDipoles(self, context)

Get the induced dipole moments of all particles.

Parameters:

context (Context) – the Context for which to get the induced dipoles

Returns:

dipoles – the induced dipole moment of particle i is stored into the i’th element

Return type:

Sequence[Vec3]

getTotalDipoles(self, context)

Get the total dipole moments (fixed plus induced) of all particles.

Parameters:

context (Context) – the Context for which to get the total dipoles

Returns:

dipoles – the total dipole moment of particle i is stored into the i’th element

Return type:

Sequence[Vec3]

getElectrostaticPotential(self, inputGrid, context)

Get the electrostatic potential.

Parameters:
  • inputGrid (Sequence[Vec3]) – input grid points over which the potential is to be evaluated

  • context (Context) – context

Returns:

outputElectrostaticPotential – output potential

Return type:

Sequence[float]

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.

Parameters:

context (Context) – context

Returns:

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)

Return type:

Sequence[float]

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.

usesPeriodicBoundaryConditions(self) bool

Returns whether or not this force makes use of periodic boundary conditions.

Returns:

true if nonbondedMethod uses PBC and false otherwise

Return type:

bool

getForceGroup(self) int

Get the force group this Force belongs to.

getName(self) str

Get the name of this Force. This is an arbitrary, user modifiable identifier. By default it equals the class name, but you can change it to anything useful.

setForceGroup(self, group)

Set the force group this Force belongs to.

Parameters:

group (int) – the group index. Legal values are between 0 and 31 (inclusive).

setName(self, name)

Set the name of this Force. This is an arbitrary, user modifiable identifier. By default it equals the class name, but you can change it to anything useful.