Class GeneralizedKirkwoodForce

java.lang.Object
ffx.openmm.Force
ffx.openmm.amoeba.GeneralizedKirkwoodForce
Direct Known Subclasses:
AmoebaGeneralizedKirkwoodForce

public class GeneralizedKirkwoodForce extends Force
This class implements an implicit solvation force using the Amoeba Generalized Kirkwood model.

To use this class, create a GeneralizedKirkwoodForce object, then call addParticle() once for each particle in the System to define its parameters. The number of particles for which you define parameters must be equal to the number of particles in the System, or else an exception will be thrown when you try to create a Context. After a particle 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().

The force supports both 3-parameter and 5-parameter particle definitions, where the extended form includes additional descreening and neck correction parameters for enhanced accuracy in specific molecular environments.

  • Constructor Details

    • GeneralizedKirkwoodForce

      public GeneralizedKirkwoodForce()
  • Method Details

    • addParticle

      public int addParticle(double charge, double radius, double scalingFactor)
      Add the 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.

      This method is provided for backwards compatibility. Compared to the alternative five parameter addParticle method, the descreenRadius parameter is set to base radius value and the neckFactor is set to zero (no neck descreening).

      Parameters:
      charge - the charge of the particle, measured in units of the proton charge
      radius - the atomic radius of the particle, measured in nm
      scalingFactor - the scaling factor for the particle
      Returns:
      the index of the particle that was added
    • addParticle

      public int addParticle(double charge, double radius, double scalingFactor, double descreenRadius, double neckFactor)
      Add the 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.

      For generalized Born / generalized Kirkwood methods, the radius of each atom has two roles. The first is to define the base radius of the atom when computing its effective radius. This base radius is usually parameterized against solvation free energy differences. The second role is to describe how much continuum water is displaced when the atom descreens water for the calculation of the Born radii of other atoms. Separation of the two roles into the "radius" and "descreenRadius" parameters gives model developers more control over these separate roles.

      For example, the fitting of base "radius" values will usually result in deviation from the force field's van der Waals definition of Rmin (or sigma). The descreenRadius can be defined separately using force field van der Waals Rmin values, which maintains consistency of atomic sizes during the HCT pairwise descreening integral. The "scalingFactor" is applied to the descreenRadius during the HCT pairwise descreening integral, while the neckFactor (if greater than zero) includes neck contributions to descreening.

      Parameters:
      charge - the charge of the particle, measured in units of the proton charge
      radius - the atomic radius of the particle, measured in nm
      scalingFactor - the scaling factor for the particle (unitless)
      descreenRadius - the atomic radius of the particle for descreening, measure in nm
      neckFactor - the scaling factor for interstitial neck descreening (unitless)
      Returns:
      the index of the particle that was added
    • destroy

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

      public double getDescreenOffset()
      Get the descreen offset.
      Returns:
      The descreen offset.
    • getDielectricOffset

      public double getDielectricOffset()
      Get the dielectric offset.
      Returns:
      The dielectric offset.
    • getIncludeCavityTerm

      public int getIncludeCavityTerm()
      Get the include cavity term.
      Returns:
      The include cavity term.
    • getNumParticles

      public int getNumParticles()
      Get the number of particles in the force.
      Returns:
      The number of particles.
    • getParticleParameters

      public void getParticleParameters(int index, com.sun.jna.ptr.DoubleByReference charge, com.sun.jna.ptr.DoubleByReference radius, com.sun.jna.ptr.DoubleByReference scalingFactor, com.sun.jna.ptr.DoubleByReference descreenRadius, com.sun.jna.ptr.DoubleByReference neckFactor)
      Get the force field parameters for a particle.
      Parameters:
      index - the index of the particle for which to get parameters
      charge - the charge of the particle, measured in units of the proton charge (output)
      radius - the atomic radius of the particle, measured in nm (output)
      scalingFactor - the scaling factor for the particle (output)
      descreenRadius - the atomic radius of the particle for descreening, measure in nm (output)
      neckFactor - the scaling factor for interstitial neck descreening (unitless) (output)
    • getProbeRadius

      public double getProbeRadius()
      Get the probe radius.
      Returns:
      The probe radius.
    • getSoluteDielectric

      public double getSoluteDielectric()
      Get the solute dielectric constant.
      Returns:
      The solute dielectric constant.
    • getSolventDielectric

      public double getSolventDielectric()
      Get the solvent dielectric constant.
      Returns:
      The solvent dielectric constant.
    • getSurfaceAreaFactor

      public double getSurfaceAreaFactor()
      Get the surface area factor.
      Returns:
      The surface area factor.
    • getTanhParameters

      public void getTanhParameters(com.sun.jna.ptr.DoubleByReference b0, com.sun.jna.ptr.DoubleByReference b1, com.sun.jna.ptr.DoubleByReference b2)
      Get Tanh function parameters b0, b1 and b2.
      Parameters:
      b0 - The first tanh parameter (output).
      b1 - The second tanh parameter (output).
      b2 - The third tanh parameter (output).
    • getTanhRescaling

      public int getTanhRescaling()
      Get the tanh rescaling.
      Returns:
      The tanh rescaling.
    • setDescreenOffset

      public void setDescreenOffset(double offset)
      Set the descreen offset.
      Parameters:
      offset - The descreen offset.
    • setDielectricOffset

      public void setDielectricOffset(double offset)
      Set the dielectric offset.
      Parameters:
      offset - The dielectric offset.
    • setIncludeCavityTerm

      public void setIncludeCavityTerm(int includeCavityTerm)
      Set the include cavity term.
      Parameters:
      includeCavityTerm - The include cavity term.
    • setParticleParameters

      public void setParticleParameters(int index, double charge, double radius, double scalingFactor, double descreenRadius, double neckFactor)
      Set the force field parameters for a particle.
      Parameters:
      index - the index of the particle for which to set parameters
      charge - the charge of the particle, measured in units of the proton charge
      radius - the atomic radius of the particle, measured in nm
      scalingFactor - the scaling factor for the particle
      descreenRadius - the atomic radius of the particle for descreening, measure in nm
      neckFactor - the scaling factor for interstitial neck descreening (unitless)
    • setProbeRadius

      public void setProbeRadius(double radius)
      Set the probe radius.
      Parameters:
      radius - The probe radius.
    • setSoluteDielectric

      public void setSoluteDielectric(double dielectric)
      Set the solute dielectric constant.
      Parameters:
      dielectric - The solute dielectric constant.
    • setSolventDielectric

      public void setSolventDielectric(double dielectric)
      Set the solvent dielectric constant.
      Parameters:
      dielectric - The solvent dielectric constant.
    • setSurfaceAreaFactor

      public void setSurfaceAreaFactor(double surfaceAreaFactor)
      Set the surface area factor.
      Parameters:
      surfaceAreaFactor - The surface area factor.
    • setTanhParameters

      public void setTanhParameters(double beta0, double beta1, double beta2)
      Set the tanh parameters.
      Parameters:
      beta0 - The tanh parameter beta0.
      beta1 - The tanh parameter beta1.
      beta2 - The tanh parameter beta2.
    • setTanhRescaling

      public void setTanhRescaling(int tanhRescale)
      Set the tanh rescaling.
      Parameters:
      tanhRescale - The tanh rescaling.
    • 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.