OpenMM
|
This class implements the Amoeba multipole interaction. More...
#include <AmoebaMultipoleForce.h>
Public Types | |
enum | NonbondedMethod { NoCutoff = 0, PME = 1 } |
enum | PolarizationType { Mutual = 0, Direct = 1 } |
enum | MultipoleAxisTypes { ZThenX = 0, Bisector = 1, ZBisect = 2, ThreeFold = 3, ZOnly = 4, NoAxisType = 5, LastAxisTypeIndex = 6 } |
enum | CovalentType { Covalent12 = 0, Covalent13 = 1, Covalent14 = 2, Covalent15 = 3, PolarizationCovalent11 = 4, PolarizationCovalent12 = 5, PolarizationCovalent13 = 6, PolarizationCovalent14 = 7, CovalentEnd = 8 } |
Public Member Functions | |
AmoebaMultipoleForce () | |
Create an AmoebaMultipoleForce. More... | |
int | getNumMultipoles () const |
Get the number of particles in the potential function. More... | |
NonbondedMethod | getNonbondedMethod () const |
Get the method used for handling long-range nonbonded interactions. More... | |
void | setNonbondedMethod (NonbondedMethod method) |
Set the method used for handling long-range nonbonded interactions. More... | |
PolarizationType | getPolarizationType () const |
Get polarization type. More... | |
void | setPolarizationType (PolarizationType type) |
Set the polarization type. More... | |
double | getCutoffDistance () const |
Get the cutoff distance (in nm) being used for nonbonded interactions. More... | |
void | setCutoffDistance (double distance) |
Set the cutoff distance (in nm) being used for nonbonded interactions. More... | |
double | getAEwald () const |
Get the Ewald alpha parameter. More... | |
void | setAEwald (double aewald) |
Set the Ewald alpha parameter. More... | |
int | getPmeBSplineOrder () const |
Get the B-spline order to use for PME charge spreading. More... | |
void | getPmeGridDimensions (std::vector< int > &gridDimension) const |
Get the PME grid dimensions. More... | |
void | setPmeGridDimensions (const std::vector< int > &gridDimension) |
Set the PME grid dimensions. More... | |
int | addMultipole (double charge, const std::vector< double > &molecularDipole, const std::vector< double > &molecularQuadrupole, int axisType, int multipoleAtomZ, int multipoleAtomX, int multipoleAtomY, double thole, double dampingFactor, double polarity) |
Add multipole-related info for a particle. More... | |
void | getMultipoleParameters (int index, double &charge, std::vector< double > &molecularDipole, std::vector< double > &molecularQuadrupole, int &axisType, int &multipoleAtomZ, int &multipoleAtomX, int &multipoleAtomY, double &thole, double &dampingFactor, double &polarity) const |
Get the multipole parameters for a particle. More... | |
void | setMultipoleParameters (int index, double charge, const std::vector< double > &molecularDipole, const std::vector< double > &molecularQuadrupole, int axisType, int multipoleAtomZ, int multipoleAtomX, int multipoleAtomY, double thole, double dampingFactor, double polarity) |
Set the multipole parameters for a particle. More... | |
void | setCovalentMap (int index, CovalentType typeId, const std::vector< int > &covalentAtoms) |
Set the CovalentMap for an atom. More... | |
void | getCovalentMap (int index, CovalentType typeId, std::vector< int > &covalentAtoms) const |
Get the CovalentMap for an atom. More... | |
void | getCovalentMaps (int index, std::vector< std::vector< int > > &covalentLists) const |
Get the CovalentMap for an atom. More... | |
int | getMutualInducedMaxIterations (void) const |
Get the max number of iterations to be used in calculating the mutual induced dipoles. More... | |
void | setMutualInducedMaxIterations (int inputMutualInducedMaxIterations) |
Set the max number of iterations to be used in calculating the mutual induced dipoles. More... | |
double | getMutualInducedTargetEpsilon (void) const |
Get the target epsilon to be used to test for convergence of iterative method used in calculating the mutual induced dipoles. More... | |
void | setMutualInducedTargetEpsilon (double inputMutualInducedTargetEpsilon) |
Set the target epsilon to be used to test for convergence of iterative method used in calculating the mutual induced dipoles. More... | |
double | getEwaldErrorTolerance () const |
Get the error tolerance for Ewald summation. More... | |
void | setEwaldErrorTolerance (double tol) |
Get the error tolerance for Ewald summation. More... | |
void | getInducedDipoles (Context &context, std::vector< Vec3 > &dipoles) |
Get the induced dipole moments of all particles. More... | |
void | getElectrostaticPotential (const std::vector< Vec3 > &inputGrid, Context &context, std::vector< double > &outputElectrostaticPotential) |
Get the electrostatic potential. More... | |
void | getSystemMultipoleMoments (Context &context, std::vector< double > &outputMultipoleMoments) |
Get the system multipole moments. More... | |
void | updateParametersInContext (Context &context) |
Update the multipole parameters in a Context to match those stored in this Force object. More... | |
bool | usesPeriodicBoundaryConditions () const |
Returns whether or not this force makes use of periodic boundary conditions. More... | |
![]() | |
Force () | |
virtual | ~Force () |
int | getForceGroup () const |
Get the force group this Force belongs to. More... | |
void | setForceGroup (int group) |
Set the force group this Force belongs to. More... | |
Protected Member Functions | |
ForceImpl * | createImpl () const |
When a Context is created, it invokes this method on each Force in the System. More... | |
![]() | |
ForceImpl & | getImplInContext (Context &context) |
Get the ForceImpl corresponding to this Force in a Context. More... | |
ContextImpl & | getContextImpl (Context &context) |
Get the ContextImpl corresponding to a Context. More... | |
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().
enum CovalentType |
enum MultipoleAxisTypes |
enum NonbondedMethod |
enum PolarizationType |
Create an AmoebaMultipoleForce.
int addMultipole | ( | double | charge, |
const std::vector< double > & | molecularDipole, | ||
const std::vector< double > & | molecularQuadrupole, | ||
int | axisType, | ||
int | multipoleAtomZ, | ||
int | multipoleAtomX, | ||
int | multipoleAtomY, | ||
double | thole, | ||
double | dampingFactor, | ||
double | polarity | ||
) |
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 |
|
protectedvirtual |
double getAEwald | ( | ) | const |
Get the Ewald alpha parameter.
If this is 0 (the default), a value is chosen automatically based on the Ewald error tolerance.
void getCovalentMap | ( | int | index, |
CovalentType | typeId, | ||
std::vector< int > & | covalentAtoms | ||
) | const |
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 |
void getCovalentMaps | ( | int | index, |
std::vector< std::vector< int > > & | covalentLists | ||
) | const |
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 |
double getCutoffDistance | ( | ) | const |
Get the cutoff distance (in nm) being used for nonbonded interactions.
If the NonbondedMethod in use is NoCutoff, this value will have no effect.
void getElectrostaticPotential | ( | const std::vector< Vec3 > & | inputGrid, |
Context & | context, | ||
std::vector< double > & | outputElectrostaticPotential | ||
) |
Get the electrostatic potential.
inputGrid | input grid points over which the potential is to be evaluated |
context | context |
outputElectrostaticPotential | output potential |
double getEwaldErrorTolerance | ( | ) | const |
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.
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 |
void getMultipoleParameters | ( | int | index, |
double & | charge, | ||
std::vector< double > & | molecularDipole, | ||
std::vector< double > & | molecularQuadrupole, | ||
int & | axisType, | ||
int & | multipoleAtomZ, | ||
int & | multipoleAtomX, | ||
int & | multipoleAtomY, | ||
double & | thole, | ||
double & | dampingFactor, | ||
double & | polarity | ||
) | const |
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 |
int getMutualInducedMaxIterations | ( | void | ) | const |
Get the max number of iterations to be used in calculating the mutual induced dipoles.
double getMutualInducedTargetEpsilon | ( | void | ) | const |
Get the target epsilon to be used to test for convergence of iterative method used in calculating the mutual induced dipoles.
NonbondedMethod getNonbondedMethod | ( | ) | const |
Get the method used for handling long-range nonbonded interactions.
|
inline |
Get the number of particles in the potential function.
int getPmeBSplineOrder | ( | ) | const |
Get the B-spline order to use for PME charge spreading.
void getPmeGridDimensions | ( | std::vector< int > & | gridDimension | ) | const |
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.
PolarizationType getPolarizationType | ( | ) | const |
Get polarization type.
void getSystemMultipoleMoments | ( | Context & | context, |
std::vector< double > & | outputMultipoleMoments | ||
) |
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) |
void setAEwald | ( | 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 |
void setCovalentMap | ( | int | index, |
CovalentType | typeId, | ||
const std::vector< int > & | 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 |
void setCutoffDistance | ( | 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 |
void setEwaldErrorTolerance | ( | 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.
void setMultipoleParameters | ( | int | index, |
double | charge, | ||
const std::vector< double > & | molecularDipole, | ||
const std::vector< double > & | 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 |
void setMutualInducedMaxIterations | ( | int | inputMutualInducedMaxIterations | ) |
Set the max number of iterations to be used in calculating the mutual induced dipoles.
max | number of iterations |
void setMutualInducedTargetEpsilon | ( | 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 |
void setNonbondedMethod | ( | NonbondedMethod | method | ) |
Set the method used for handling long-range nonbonded interactions.
void setPmeGridDimensions | ( | const std::vector< int > & | 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 |
void setPolarizationType | ( | PolarizationType | type | ) |
Set the polarization type.
void updateParametersInContext | ( | 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 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.
|
inlinevirtual |
Returns whether or not this force makes use of periodic boundary conditions.
Reimplemented from Force.
References AmoebaMultipoleForce::PME.