Class MultipoleForce

java.lang.Object
ffx.openmm.Force
ffx.openmm.amoeba.MultipoleForce
Direct Known Subclasses:
AmoebaMultipoleForce

public class MultipoleForce extends Force
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().

  • Constructor Details

    • MultipoleForce

      public MultipoleForce()
  • Method Details

    • addMultipole

      public int addMultipole(double charge, DoubleArray molecularDipole, DoubleArray molecularQuadrupole, int axisType, int multipoleAtomZ, int multipoleAtomX, int multipoleAtomY, double thole, double dampingFactor, double polarity)
      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
      multipoleAtomZ - index of first atom used in constructing lab to molecular frames
      multipoleAtomX - index of second atom used in constructing lab to molecular frames
      multipoleAtomY - index of second atom used in constructing lab to molecular frames
      thole - Thole parameter
      dampingFactor - dampingFactor parameter
      polarity - polarity parameter
      Returns:
      the index of the particle that was added
    • destroy

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

      @Deprecated public double getAEwald()
      Deprecated.
      This method exists only for backward compatibility. Use getPMEParameters() instead.
      Get the Ewald alpha parameter. If this is 0 (the default), a value is chosen automatically based on the Ewald error tolerance.
      Returns:
      the Ewald alpha parameter
    • getCovalentMap

      public IntArray getCovalentMap(int i, int covalentType)
      Get the covalent map for a given atom index and covalent type.
      Parameters:
      i - The atom index.
      covalentType - The covalent type.
      Returns:
      An IntArray representing the covalent map for the specified atom index and covalent type.
    • getCovalentMaps

      public IntArray getCovalentMaps(int i)
      Get all covalent maps for a given atom index.
      Parameters:
      i - The atom index.
      Returns:
      An IntArray containing all covalent maps for the specified atom.
    • getCutoffDistance

      public double getCutoffDistance()
      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
    • getElectrostaticPotential

      public DoubleArray getElectrostaticPotential(Context context, DoubleArray points)
      Get the electrostatic potential at specified points.
      Parameters:
      context - The OpenMM context.
      points - The points at which to calculate the potential.
      Returns:
      A DoubleArray containing the electrostatic potential at each point.
    • getEwaldErrorTolerance

      public double getEwaldErrorTolerance()
      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.

    • getExtrapolationCoefficients

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

      public void getInducedDipoles(Context context, DoubleArray inducedDipoles)
      Get the induced dipoles.
      Parameters:
      context - The OpenMM context.
      inducedDipoles - The induced dipoles (output).
    • getLabFramePermanentDipoles

      public void getLabFramePermanentDipoles(Context context, DoubleArray dipoles)
      Get the lab frame permanent dipoles.
      Parameters:
      context - The OpenMM context.
      dipoles - The lab frame permanent dipoles (output).
    • getMultipoleParameters

      public void getMultipoleParameters(int index, com.sun.jna.ptr.DoubleByReference charge, com.sun.jna.ptr.PointerByReference molecularDipole, com.sun.jna.ptr.PointerByReference molecularQuadrupole, com.sun.jna.ptr.IntByReference axisType, com.sun.jna.ptr.IntByReference multipoleAtomZ, com.sun.jna.ptr.IntByReference multipoleAtomX, com.sun.jna.ptr.IntByReference multipoleAtomY, com.sun.jna.ptr.DoubleByReference thole, com.sun.jna.ptr.DoubleByReference dampingFactor, com.sun.jna.ptr.DoubleByReference polarity)
      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
      multipoleAtomZ - index of first atom used in constructing lab to molecular frames
      multipoleAtomX - index of second atom used in constructing lab to molecular frames
      multipoleAtomY - index of second atom used in constructing lab to molecular frames
      thole - Thole parameter
      dampingFactor - dampingFactor parameter
      polarity - polarity parameter
    • getMutualInducedMaxIterations

      public int getMutualInducedMaxIterations()
      Get the mutual induced max iterations.
      Returns:
      The mutual induced max iterations.
    • getMutualInducedTargetEpsilon

      public double getMutualInducedTargetEpsilon()
      Get the mutual induced target epsilon.
      Returns:
      The mutual induced target epsilon.
    • getNonbondedMethod

      public int getNonbondedMethod()
      Get the nonbonded method for the multipole force.
      Returns:
      The nonbonded method.
    • getNumMultipoles

      public int getNumMultipoles()
      Get the number of multipoles.
      Returns:
      The number of multipoles.
    • 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 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
    • 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 the context.
      Parameters:
      context - The OpenMM context.
      alpha - The Ewald alpha parameter (output).
      nx - The PME grid size in x (output).
      ny - The PME grid size in y (output).
      nz - The PME grid size in z (output).
    • getPmeGridDimensions

      public IntArray getPmeGridDimensions()
      Get the PME grid dimensions.
      Returns:
      An IntArray containing the PME grid dimensions.
    • getPmeBSplineOrder

      public int getPmeBSplineOrder()
      Get the PME B-spline order.
      Returns:
      The PME B-spline order.
    • getPolarizationType

      public int getPolarizationType()
      Get the polarization type.
      Returns:
      The polarization type.
    • getSystemMultipoleMoments

      public void getSystemMultipoleMoments(Context context, DoubleArray moments)
      Get the system multipole moments.
      Parameters:
      context - The OpenMM context.
      moments - The system multipole moments (output).
    • getTotalDipoles

      public void getTotalDipoles(Context context, DoubleArray dipoles)
      Get the total dipoles.
      Parameters:
      context - The OpenMM context.
      dipoles - The total dipoles (output).
    • setAEwald

      @Deprecated public void setAEwald(double aewald)
      Deprecated.
      This method exists only for backward compatibility. Use setPMEParameters() instead.
      Set the Ewald alpha parameter. If this is 0 (the default), a value is chosen automatically based on the Ewald error tolerance.
      Parameters:
      aewald - alpha parameter
    • setCovalentMap

      public void setCovalentMap(int i, int covalentType, IntArray covalentMap)
      Set the covalent map.
      Parameters:
      i - The atom index.
      covalentType - The covalent type.
      covalentMap - The covalent map.
    • setCutoffDistance

      public void setCutoffDistance(double cutoff)
      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:
      cutoff - the cutoff distance, measured in nm
    • setEwaldErrorTolerance

      public void setEwaldErrorTolerance(double ewaldTolerance)
      Set 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.

      Parameters:
      ewaldTolerance - the error tolerance
    • setExtrapolationCoefficients

      public void setExtrapolationCoefficients(DoubleArray exptCoefficients)
      Set extrapolation coefficients.
      Parameters:
      exptCoefficients - The extrapolation coefficients.
    • setMultipoleParameters

      public void setMultipoleParameters(int index, double charge, DoubleArray molecularDipole, DoubleArray 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
      multipoleAtomZ - index of first atom used in constructing lab to molecular frames
      multipoleAtomX - index of second atom used in constructing lab to molecular frames
      multipoleAtomY - index of second atom used in constructing lab to molecular frames
      thole - thole parameter
      dampingFactor - damping factor parameter
      polarity - polarity parameter
    • setMutualInducedMaxIterations

      public void setMutualInducedMaxIterations(int iterations)
      Set the mutual induced target maximum number of iterations.
      Parameters:
      iterations - The mutual induced max iterations.
    • setMutualInducedTargetEpsilon

      public void setMutualInducedTargetEpsilon(double epsilon)
      Set the mutual induced target epsilon.
      Parameters:
      epsilon - The mutual induced target epsilon.
    • setNonbondedMethod

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

      public 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
    • setPmeGridDimensions

      public void setPmeGridDimensions(IntArray gridDimensions)
      Set the PME grid dimensions for the multipole force.
      Parameters:
      gridDimensions - The PME grid dimensions.
    • setPolarizationType

      public void setPolarizationType(int method)
      Set the polarization method.
      Parameters:
      method - The polarization method.
    • updateParametersInContext

      public void updateParametersInContext(Context context)
      Update the parameters in the context.
      Parameters:
      context - The OpenMM context.
    • 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.