Class HippoNonbondedForce

java.lang.Object
ffx.openmm.Force
ffx.openmm.amoeba.HippoNonbondedForce

public class HippoNonbondedForce extends 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.

  • Field Summary

    Fields inherited from class ffx.openmm.Force

    pointer
  • Constructor Summary

    Constructors
    Constructor
    Description
    Create a new HippoNonbondedForce.
  • Method Summary

    Modifier and Type
    Method
    Description
    int
    addException(int particle1, int particle2, double multipoleMultipoleScale, double dipoleMultipoleScale, double dipoleDipoleScale, double dispersionScale, double repulsionScale, double chargeTransferScale, int replace)
    Add an interaction to the list of exceptions that should be calculated differently from other interactions.
    int
    addParticle(double charge, com.sun.jna.ptr.PointerByReference dipole, com.sun.jna.ptr.PointerByReference 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.
    void
    Destroy the force.
    double
    Get the cutoff distance.
    void
    getDPMEParameters(com.sun.jna.ptr.DoubleByReference alpha, com.sun.jna.ptr.IntByReference nx, com.sun.jna.ptr.IntByReference ny, com.sun.jna.ptr.IntByReference nz)
    Get the DPME parameters.
    void
    Get the DPME parameters.
    void
    getDPMEParametersInContext(Context context, com.sun.jna.ptr.DoubleByReference alpha, com.sun.jna.ptr.IntByReference nx, com.sun.jna.ptr.IntByReference ny, com.sun.jna.ptr.IntByReference nz)
    Get the DPME parameters in context.
    void
    Get the DPME parameters in context.
    double
    Get the Ewald error tolerance.
    void
    getExceptionParameters(int index, com.sun.jna.ptr.IntByReference particle1, com.sun.jna.ptr.IntByReference particle2, com.sun.jna.ptr.DoubleByReference multipoleMultipoleScale, com.sun.jna.ptr.DoubleByReference dipoleMultipoleScale, com.sun.jna.ptr.DoubleByReference dipoleDipoleScale, com.sun.jna.ptr.DoubleByReference dispersionScale, com.sun.jna.ptr.DoubleByReference repulsionScale, com.sun.jna.ptr.DoubleByReference chargeTransferScale)
    Get the scale factors for an interaction that should be calculated differently from others.
    void
    getExceptionParameters(int index, IntBuffer particle1, IntBuffer particle2, DoubleBuffer multipoleMultipoleScale, DoubleBuffer dipoleMultipoleScale, DoubleBuffer dipoleDipoleScale, DoubleBuffer dispersionScale, DoubleBuffer repulsionScale, DoubleBuffer chargeTransferScale)
    Get the scale factors for an interaction that should be calculated differently from others.
    com.sun.jna.ptr.PointerByReference
    Get the extrapolation coefficients.
    void
    getInducedDipoles(Context context, com.sun.jna.ptr.PointerByReference dipoles)
    Get the induced dipoles.
    void
    getLabFramePermanentDipoles(Context context, com.sun.jna.ptr.PointerByReference dipoles)
    Get the lab frame permanent dipoles.
    int
    Get the nonbonded method.
    int
    Get the number of exceptions.
    int
    Get the number of particles.
    void
    getParticleParameters(int index, com.sun.jna.ptr.DoubleByReference charge, com.sun.jna.ptr.PointerByReference dipole, com.sun.jna.ptr.PointerByReference quadrupole, com.sun.jna.ptr.DoubleByReference coreCharge, com.sun.jna.ptr.DoubleByReference alpha, com.sun.jna.ptr.DoubleByReference epsilon, com.sun.jna.ptr.DoubleByReference damping, com.sun.jna.ptr.DoubleByReference c6, com.sun.jna.ptr.DoubleByReference pauliK, com.sun.jna.ptr.DoubleByReference pauliQ, com.sun.jna.ptr.DoubleByReference pauliAlpha, com.sun.jna.ptr.DoubleByReference polarizability, com.sun.jna.ptr.IntByReference axisType, com.sun.jna.ptr.IntByReference multipoleAtomZ, com.sun.jna.ptr.IntByReference multipoleAtomX, com.sun.jna.ptr.IntByReference multipoleAtomY)
    Get the nonbonded force parameters for a particle.
    void
    getParticleParameters(int index, DoubleBuffer charge, com.sun.jna.ptr.PointerByReference dipole, com.sun.jna.ptr.PointerByReference quadrupole, DoubleBuffer coreCharge, DoubleBuffer alpha, DoubleBuffer epsilon, DoubleBuffer damping, DoubleBuffer c6, DoubleBuffer pauliK, DoubleBuffer pauliQ, DoubleBuffer pauliAlpha, DoubleBuffer polarizability, IntBuffer axisType, IntBuffer multipoleAtomZ, IntBuffer multipoleAtomX, IntBuffer multipoleAtomY)
    Get the nonbonded force parameters for a particle.
    void
    getPMEParameters(com.sun.jna.ptr.DoubleByReference alpha, com.sun.jna.ptr.IntByReference nx, com.sun.jna.ptr.IntByReference ny, com.sun.jna.ptr.IntByReference nz)
    Get the PME parameters.
    void
    Get the PME parameters.
    void
    getPMEParametersInContext(Context context, com.sun.jna.ptr.DoubleByReference alpha, com.sun.jna.ptr.IntByReference nx, com.sun.jna.ptr.IntByReference ny, com.sun.jna.ptr.IntByReference nz)
    Get the PME parameters in context.
    void
    Get the PME parameters in context.
    double
    Get the switching distance.
    void
    setCutoffDistance(double distance)
    Set the cutoff distance.
    void
    setDPMEParameters(double alpha, int nx, int ny, int nz)
    Set the DPME parameters.
    void
    setEwaldErrorTolerance(double tolerance)
    Set the Ewald error tolerance.
    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.
    void
    setExtrapolationCoefficients(com.sun.jna.ptr.PointerByReference coefficients)
    Set the extrapolation coefficients.
    void
    setNonbondedMethod(int method)
    Set the nonbonded method.
    void
    setParticleParameters(int index, double charge, com.sun.jna.ptr.PointerByReference dipole, com.sun.jna.ptr.PointerByReference 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.
    void
    setPMEParameters(double alpha, int nx, int ny, int nz)
    Set the PME parameters.
    void
    setSwitchingDistance(double distance)
    Set the switching distance.
    void
    Update the parameters in a Context to match those stored in this Force object.
    boolean
    Check if the force uses periodic boundary conditions.

    Methods inherited from class ffx.openmm.Force

    getForceGroup, getForceIndex, getName, getPointer, setForceGroup, setForceIndex, setName

    Methods inherited from class java.lang.Object

    clone, equals, finalize, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
  • Constructor Details

    • HippoNonbondedForce

      public HippoNonbondedForce()
      Create a new HippoNonbondedForce.
  • Method Details

    • addException

      public int addException(int particle1, int particle2, double multipoleMultipoleScale, double dipoleMultipoleScale, double dipoleDipoleScale, double dispersionScale, double repulsionScale, double chargeTransferScale, int replace)
      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
    • addParticle

      public int addParticle(double charge, com.sun.jna.ptr.PointerByReference dipole, com.sun.jna.ptr.PointerByReference 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
    • destroy

      public void destroy()
      Destroy the force.
      Specified by:
      destroy in class Force
    • getCutoffDistance

      public double getCutoffDistance()
      Get the cutoff distance.
      Returns:
      The cutoff distance, measured in nm.
    • getDPMEParameters

      public void getDPMEParameters(com.sun.jna.ptr.DoubleByReference alpha, com.sun.jna.ptr.IntByReference nx, com.sun.jna.ptr.IntByReference ny, com.sun.jna.ptr.IntByReference nz)
      Get the DPME parameters.
      Parameters:
      alpha - The Ewald alpha parameter (output).
      nx - The number of grid points along the x axis (output).
      ny - The number of grid points along the y axis (output).
      nz - The number of grid points along the z axis (output).
    • getDPMEParameters

      public void getDPMEParameters(DoubleBuffer alpha, IntBuffer nx, IntBuffer ny, IntBuffer nz)
      Get the DPME parameters.
      Parameters:
      alpha - The Ewald alpha parameter (output).
      nx - The number of grid points along the x axis (output).
      ny - The number of grid points along the y axis (output).
      nz - The number of grid points along the z axis (output).
    • getDPMEParametersInContext

      public void getDPMEParametersInContext(Context context, com.sun.jna.ptr.DoubleByReference alpha, com.sun.jna.ptr.IntByReference nx, com.sun.jna.ptr.IntByReference ny, com.sun.jna.ptr.IntByReference nz)
      Get the DPME parameters in context.
      Parameters:
      context - The context.
      alpha - The Ewald alpha parameter (output).
      nx - The number of grid points along the x axis (output).
      ny - The number of grid points along the y axis (output).
      nz - The number of grid points along the z axis (output).
    • getDPMEParametersInContext

      public void getDPMEParametersInContext(Context context, DoubleBuffer alpha, IntBuffer nx, IntBuffer ny, IntBuffer nz)
      Get the DPME parameters in context.
      Parameters:
      context - The context.
      alpha - The Ewald alpha parameter (output).
      nx - The number of grid points along the x axis (output).
      ny - The number of grid points along the y axis (output).
      nz - The number of grid points along the z axis (output).
    • getEwaldErrorTolerance

      public double getEwaldErrorTolerance()
      Get the Ewald error tolerance.
      Returns:
      The Ewald error tolerance.
    • getExceptionParameters

      public void getExceptionParameters(int index, com.sun.jna.ptr.IntByReference particle1, com.sun.jna.ptr.IntByReference particle2, com.sun.jna.ptr.DoubleByReference multipoleMultipoleScale, com.sun.jna.ptr.DoubleByReference dipoleMultipoleScale, com.sun.jna.ptr.DoubleByReference dipoleDipoleScale, com.sun.jna.ptr.DoubleByReference dispersionScale, com.sun.jna.ptr.DoubleByReference repulsionScale, com.sun.jna.ptr.DoubleByReference chargeTransferScale)
      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
    • getExceptionParameters

      public void getExceptionParameters(int index, IntBuffer particle1, IntBuffer particle2, DoubleBuffer multipoleMultipoleScale, DoubleBuffer dipoleMultipoleScale, DoubleBuffer dipoleDipoleScale, DoubleBuffer dispersionScale, DoubleBuffer repulsionScale, DoubleBuffer chargeTransferScale)
      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
    • getExtrapolationCoefficients

      public com.sun.jna.ptr.PointerByReference getExtrapolationCoefficients()
      Get the extrapolation coefficients.
      Returns:
      The extrapolation coefficients.
    • getInducedDipoles

      public void getInducedDipoles(Context context, com.sun.jna.ptr.PointerByReference dipoles)
      Get the induced dipoles.
      Parameters:
      context - The context.
      dipoles - The induced dipoles (output).
    • getLabFramePermanentDipoles

      public void getLabFramePermanentDipoles(Context context, com.sun.jna.ptr.PointerByReference dipoles)
      Get the lab frame permanent dipoles.
      Parameters:
      context - The context.
      dipoles - The lab frame permanent dipoles (output).
    • getNonbondedMethod

      public int getNonbondedMethod()
      Get the nonbonded method.
      Returns:
      The nonbonded method.
    • getNumExceptions

      public int getNumExceptions()
      Get the number of exceptions.
      Returns:
      The number of exceptions.
    • getNumParticles

      public int getNumParticles()
      Get the number of particles.
      Returns:
      The number of particles.
    • getPMEParameters

      public void getPMEParameters(com.sun.jna.ptr.DoubleByReference alpha, com.sun.jna.ptr.IntByReference nx, com.sun.jna.ptr.IntByReference ny, com.sun.jna.ptr.IntByReference nz)
      Get the PME parameters.
      Parameters:
      alpha - The Ewald alpha parameter (output).
      nx - The number of grid points along the x axis (output).
      ny - The number of grid points along the y axis (output).
      nz - The number of grid points along the z axis (output).
    • getPMEParameters

      public void getPMEParameters(DoubleBuffer alpha, IntBuffer nx, IntBuffer ny, IntBuffer nz)
      Get the PME parameters.
      Parameters:
      alpha - The Ewald alpha parameter (output).
      nx - The number of grid points along the x axis (output).
      ny - The number of grid points along the y axis (output).
      nz - The number of grid points along the z axis (output).
    • getPMEParametersInContext

      public void getPMEParametersInContext(Context context, com.sun.jna.ptr.DoubleByReference alpha, com.sun.jna.ptr.IntByReference nx, com.sun.jna.ptr.IntByReference ny, com.sun.jna.ptr.IntByReference nz)
      Get the PME parameters in context.
      Parameters:
      context - The context.
      alpha - The Ewald alpha parameter (output).
      nx - The number of grid points along the x axis (output).
      ny - The number of grid points along the y axis (output).
      nz - The number of grid points along the z axis (output).
    • getPMEParametersInContext

      public void getPMEParametersInContext(Context context, DoubleBuffer alpha, IntBuffer nx, IntBuffer ny, IntBuffer nz)
      Get the PME parameters in context.
      Parameters:
      context - The context.
      alpha - The Ewald alpha parameter (output).
      nx - The number of grid points along the x axis (output).
      ny - The number of grid points along the y axis (output).
      nz - The number of grid points along the z axis (output).
    • getParticleParameters

      public void getParticleParameters(int index, com.sun.jna.ptr.DoubleByReference charge, com.sun.jna.ptr.PointerByReference dipole, com.sun.jna.ptr.PointerByReference quadrupole, com.sun.jna.ptr.DoubleByReference coreCharge, com.sun.jna.ptr.DoubleByReference alpha, com.sun.jna.ptr.DoubleByReference epsilon, com.sun.jna.ptr.DoubleByReference damping, com.sun.jna.ptr.DoubleByReference c6, com.sun.jna.ptr.DoubleByReference pauliK, com.sun.jna.ptr.DoubleByReference pauliQ, com.sun.jna.ptr.DoubleByReference pauliAlpha, com.sun.jna.ptr.DoubleByReference polarizability, com.sun.jna.ptr.IntByReference axisType, com.sun.jna.ptr.IntByReference multipoleAtomZ, com.sun.jna.ptr.IntByReference multipoleAtomX, com.sun.jna.ptr.IntByReference multipoleAtomY)
      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
    • getParticleParameters

      public void getParticleParameters(int index, DoubleBuffer charge, com.sun.jna.ptr.PointerByReference dipole, com.sun.jna.ptr.PointerByReference quadrupole, DoubleBuffer coreCharge, DoubleBuffer alpha, DoubleBuffer epsilon, DoubleBuffer damping, DoubleBuffer c6, DoubleBuffer pauliK, DoubleBuffer pauliQ, DoubleBuffer pauliAlpha, DoubleBuffer polarizability, IntBuffer axisType, IntBuffer multipoleAtomZ, IntBuffer multipoleAtomX, IntBuffer multipoleAtomY)
      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
    • getSwitchingDistance

      public double getSwitchingDistance()
      Get the switching distance.
      Returns:
      The switching distance, measured in nm.
    • setCutoffDistance

      public void setCutoffDistance(double distance)
      Set the cutoff distance.
      Parameters:
      distance - The cutoff distance, measured in nm.
    • setDPMEParameters

      public void setDPMEParameters(double alpha, int nx, int ny, int nz)
      Set the DPME parameters.
      Parameters:
      alpha - The Ewald alpha 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.
    • setEwaldErrorTolerance

      public void setEwaldErrorTolerance(double tolerance)
      Set the Ewald error tolerance.
      Parameters:
      tolerance - The Ewald error tolerance.
    • setExceptionParameters

      public 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.
      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
    • setExtrapolationCoefficients

      public void setExtrapolationCoefficients(com.sun.jna.ptr.PointerByReference coefficients)
      Set the extrapolation coefficients.
      Parameters:
      coefficients - The extrapolation coefficients.
    • setNonbondedMethod

      public void setNonbondedMethod(int method)
      Set the nonbonded method.
      Parameters:
      method - The nonbonded method.
    • setPMEParameters

      public void setPMEParameters(double alpha, int nx, int ny, int nz)
      Set the PME parameters.
      Parameters:
      alpha - The Ewald alpha 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.
    • setParticleParameters

      public void setParticleParameters(int index, double charge, com.sun.jna.ptr.PointerByReference dipole, com.sun.jna.ptr.PointerByReference 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
    • setSwitchingDistance

      public void setSwitchingDistance(double distance)
      Set the switching distance.
      Parameters:
      distance - The switching distance, measured in nm.
    • updateParametersInContext

      public void updateParametersInContext(Context context)
      Update the parameters in a Context to match those stored in this Force object.
      Parameters:
      context - The Context in which to update the parameters.
    • usesPeriodicBoundaryConditions

      public boolean usesPeriodicBoundaryConditions()
      Check if the force uses periodic boundary conditions.
      Overrides:
      usesPeriodicBoundaryConditions in class Force
      Returns:
      True if the force uses periodic boundary conditions.