HippoNonbondedForce

class HippoNonbondedForce : public OpenMM::Force

This class implements all nonbonded interactions in the HIPPO force field: electrostatics, induction, charge transfer, dispersion, and repulsion. Although some of these are conceptually distinct, they share parameters in common and are most efficiently computed together. For example, the same multipole definitions are used for both electrostatics and Pauli repulsion. Therefore, all of them are computed by a single Force object.

To use it, create a HippoNonbondedForce object, then call addParticle() once for each particle. After an entry has been added, you can modify its force field parameters by calling setParticleParameters(). This will have no effect on Contexts that already exist unless you call updateParametersInContext().

You also can specify “exceptions”, particular pairs of particles whose interactions should be reduced or completely omitted. Call addException() to define exceptions.

Public Functions

HippoNonbondedForce()

Create a HippoNonbondedForce.

inline int getNumParticles() const

Get the number of particles in the potential function.

inline int getNumExceptions() const

Get the number of exceptions.

NonbondedMethod getNonbondedMethod() const

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

void setNonbondedMethod(NonbondedMethod method)

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

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.

Returns

the cutoff distance, measured in nm

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.

Parameters

distance – the cutoff distance, measured in nm

double getSwitchingDistance() const

Get the distance at which the switching function begins to reduce the repulsion and charge transfer interactions. This must be less than the cutoff distance.

void setSwitchingDistance(double distance)

Set the distance at which the switching function begins to reduce the repulsion and charge transfer interactions. This must be less than the cutoff distance.

const std::vector<double> &getExtrapolationCoefficients() const

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

void setExtrapolationCoefficients(const std::vector<double> &coefficients)

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

Parameters

coefficients – a vector whose mth entry specifies the coefficient for mu_m. The length of this vector determines how many iterations are performed.

void getPMEParameters(double &alpha, int &nx, int &ny, int &nz) const

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.

Parameters
  • alpha[out] the separation parameter

  • nx[out] the number of grid points along the X axis

  • ny[out] the number of grid points along the Y axis

  • nz[out] the number of grid points along the Z axis

void getDPMEParameters(double &alpha, int &nx, int &ny, int &nz) const

Get the parameters to use for dispersion 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[out] the separation parameter

  • nx[out] the number of dispersion grid points along the X axis

  • ny[out] the number of dispersion grid points along the Y axis

  • nz[out] the number of dispersion grid points along the Z axis

void setPMEParameters(double alpha, int nx, int ny, int 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 – the separation parameter

  • nx – the number of grid points along the X axis

  • ny – the number of grid points along the Y axis

  • nz – the number of grid points along the Z axis

void setDPMEParameters(double alpha, int nx, int ny, int nz)

Set the parameters to use for dispersion 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 – the separation parameter

  • nx – the number of grid points along the X axis

  • ny – the number of grid points along the Y axis

  • nz – the number of grid points along the Z axis

void getPMEParametersInContext(const Context &context, double &alpha, int &nx, int &ny, int &nz) const

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 – the Context for which to get the parameters

  • alpha[out] the separation parameter

  • nx[out] the number of grid points along the X axis

  • ny[out] the number of grid points along the Y axis

  • nz[out] the number of grid points along the Z axis

void getDPMEParametersInContext(const Context &context, double &alpha, int &nx, int &ny, int &nz) const

Get the parameters being used for dispersion 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 setPMEParameters(), or the standard values calculated based on the Ewald error tolerance. See the manual for details.

Parameters
  • context – the Context for which to get the parameters

  • alpha[out] the separation parameter

  • nx[out] the number of grid points along the X axis

  • ny[out] the number of grid points along the Y axis

  • nz[out] the number of grid points along the Z axis

int addParticle(double charge, const std::vector<double> &dipole, const std::vector<double> &quadrupole, double coreCharge, double alpha, double epsilon, double damping, double c6, double pauliK, double pauliQ, double pauliAlpha, double polarizability, int axisType, int multipoleAtomZ, int multipoleAtomX, int multipoleAtomY)

Add the nonbonded force parameters for a particle. This should be called once for each particle in the System. When it is called for the i’th time, it specifies the parameters for the i’th particle.

Parameters
  • charge – the particle’s charge

  • dipole – the particle’s molecular dipole (vector of size 3)

  • quadrupole – the particle’s molecular quadrupole (vector of size 9)

  • coreCharge – the charge of the atomic core

  • alpha – controls the width of the particle’s electron density

  • epsilon – sets the magnitude of charge transfer

  • damping – sets the length scale for charge transfer

  • c6 – the coefficient of the dispersion interaction

  • pauliK – the coefficient of the Pauli repulsion interaction

  • pauliQ – the charge used in computing the Pauli repulsion interaction

  • pauliAlpha – the width of the particle’s electron density for computing the Pauli repulsion interaction

  • polarizability – atomic polarizability

  • axisType – the particle’s axis type

  • multipoleAtomZ – index of first atom used in defining the local coordinate system for multipoles

  • multipoleAtomX – index of second atom used in defining the local coordinate system for multipoles

  • multipoleAtomY – index of third atom used in defining the local coordinate system for multipoles

Returns

the index of the particle that was added

void getParticleParameters(int index, double &charge, std::vector<double> &dipole, std::vector<double> &quadrupole, double &coreCharge, double &alpha, double &epsilon, double &damping, double &c6, double &pauliK, double &pauliQ, double &pauliAlpha, double &polarizability, int &axisType, int &multipoleAtomZ, int &multipoleAtomX, int &multipoleAtomY) const

Get the nonbonded force parameters for a particle.

Parameters
  • index – the index of the particle for which to get parameters

  • charge – the particle’s charge

  • dipole – the particle’s molecular dipole (vector of size 3)

  • quadrupole – the particle’s molecular quadrupole (vector of size 9)

  • coreCharge – the charge of the atomic core

  • alpha – controls the width of the particle’s electron density

  • epsilon – sets the magnitude of charge transfer

  • damping – sets the length scale for charge transfer

  • c6 – the coefficient of the dispersion interaction

  • pauliK – the coefficient of the Pauli repulsion interaction

  • pauliQ – the charge used in computing the Pauli repulsion interaction

  • pauliAlpha – the width of the particle’s electron density for computing the Pauli repulsion interaction

  • polarizability – atomic polarizability

  • axisType – the particle’s axis type

  • multipoleAtomZ – index of first atom used in defining the local coordinate system for multipoles

  • multipoleAtomX – index of second atom used in defining the local coordinate system for multipoles

  • multipoleAtomY – index of third atom used in defining the local coordinate system for multipoles

void setParticleParameters(int index, double charge, const std::vector<double> &dipole, const std::vector<double> &quadrupole, double coreCharge, double alpha, double epsilon, double damping, double c6, double pauliK, double pauliQ, double pauliAlpha, double polarizability, int axisType, int multipoleAtomZ, int multipoleAtomX, int multipoleAtomY)

Set the nonbonded force parameters for a particle.

Parameters
  • index – the index of the particle for which to set parameters

  • charge – the particle’s charge

  • dipole – the particle’s molecular dipole (vector of size 3)

  • quadrupole – the particle’s molecular quadrupole (vector of size 9)

  • coreCharge – the charge of the atomic core

  • alpha – controls the width of the particle’s electron density

  • epsilon – sets the magnitude of charge transfer

  • damping – sets the length scale for charge transfer

  • c6 – the coefficient of the dispersion interaction

  • pauliK – the coefficient of the Pauli repulsion interaction

  • pauliQ – the charge used in computing the Pauli repulsion interaction

  • pauliAlpha – the width of the particle’s electron density for computing the Pauli repulsion interaction

  • polarizability – atomic polarizability

  • axisType – the particle’s axis type

  • multipoleAtomZ – index of first atom used in defining the local coordinate system for multipoles

  • multipoleAtomX – index of second atom used in defining the local coordinate system for multipoles

  • multipoleAtomY – index of third atom used in defining the local coordinate system for multipoles

int addException(int particle1, int particle2, double multipoleMultipoleScale, double dipoleMultipoleScale, double dipoleDipoleScale, double dispersionScale, double repulsionScale, double chargeTransferScale, bool replace = false)

Add an interaction to the list of exceptions that should be calculated differently from other interactions. If all scale factors are set to 0, this will cause the interaction to be completely omitted from force and energy calculations.

Parameters
  • particle1 – the index of the first particle involved in the interaction

  • particle2 – the index of the second particle involved in the interaction

  • multipoleMultipoleScale – the factor by which to scale the Coulomb interaction between fixed multipoles

  • dipoleMultipoleScale – the factor by which to scale the Coulomb interaction between an induced dipole and a fixed multipole

  • dipoleDipoleScale – the factor by which to scale the Coulomb interaction between induced dipoles

  • dispersionScale – the factor by which to scale the dispersion interaction

  • repulsionScale – the factor by which to scale the Pauli repulsion

  • chargeTransferScale – the factor by which to scale the charge transfer interaction

  • replace – determines the behavior if there is already an exception for the same two particles. If true, the existing one is replaced. If false, an exception is thrown.

Returns

the index of the exception that was added

void getExceptionParameters(int index, int &particle1, int &particle2, double &multipoleMultipoleScale, double &dipoleMultipoleScale, double &dipoleDipoleScale, double &dispersionScale, double &repulsionScale, double &chargeTransferScale) const

Get the scale factors for an interaction that should be calculated differently from others.

Parameters
  • index – the index of the interaction for which to get parameters

  • particle1 – the index of the first particle involved in the interaction

  • particle2 – the index of the second particle involved in the interaction

  • multipoleMultipoleScale – the factor by which to scale the Coulomb interaction between fixed multipoles

  • dipoleMultipoleScale – the factor by which to scale the Coulomb interaction between an induced dipole and a fixed multipole

  • dipoleDipoleScale – the factor by which to scale the Coulomb interaction between induced dipoles

  • dispersionScale – the factor by which to scale the dispersion interaction

  • repulsionScale – the factor by which to scale the Pauli repulsion

  • chargeTransferScale – the factor by which to scale the charge transfer interaction

void setExceptionParameters(int index, int particle1, int particle2, double multipoleMultipoleScale, double dipoleMultipoleScale, double dipoleDipoleScale, double dispersionScale, double repulsionScale, double chargeTransferScale)

Set the scale factors for an interaction that should be calculated differently from others. If all scale factors are set to 0, this will cause the interaction to be completely omitted from force and energy calculations.

Parameters
  • index – the index of the interaction for which to set parameters

  • particle1 – the index of the first particle involved in the interaction

  • particle2 – the index of the second particle involved in the interaction

  • multipoleMultipoleScale – the factor by which to scale the Coulomb interaction between fixed multipoles

  • dipoleMultipoleScale – the factor by which to scale the Coulomb interaction between an induced dipole and a fixed multipole

  • dipoleDipoleScale – the factor by which to scale the Coulomb interaction between induced dipoles

  • dispersionScale – the factor by which to scale the dispersion interaction

  • repulsionScale – the factor by which to scale the Pauli repulsion

  • chargeTransferScale – the factor by which to scale the charge transfer interaction

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.

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 getLabFramePermanentDipoles(Context &context, std::vector<Vec3> &dipoles)

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

Parameters
  • context – the Context for which to get the fixed dipoles

  • dipoles[out] the fixed dipole moment of particle i is stored into the i’th element

void getInducedDipoles(Context &context, std::vector<Vec3> &dipoles)

Get the induced dipole moments of all particles.

Parameters
  • context – the Context for which to get the induced dipoles

  • dipoles[out] the induced dipole moment of particle i is stored into the i’th element

void updateParametersInContext(Context &context)

Update the particle and exception 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 setParticleParameters() 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 particles and exceptions. 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, only the scale factors for an exception can be changed; the pair of particles involved in the exception cannot change. Finally, this method cannot be used to add new particles or exceptions, only to change the parameters of existing ones.

inline virtual bool usesPeriodicBoundaryConditions() const

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

Returns

true if nonbondedMethod uses PBC and false otherwise