Package ffx.openmm

Class GBSAOBCForce

java.lang.Object
ffx.openmm.Force
ffx.openmm.GBSAOBCForce

public class GBSAOBCForce extends Force
This class implements an implicit solvation force using the GBSA-OBC model.

To use this class, create a GBSAOBCForce object, then call addParticle() once for each particle in the System to define its parameters. The number of particles for which you define GBSA parameters must be exactly 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().

When using this Force, the System should also include a NonbondedForce, and both objects must specify identical charges for all particles. Otherwise, the results will not be correct. Furthermore, if the nonbonded method is set to CutoffNonPeriodic or CutoffPeriodic, you should call setReactionFieldDielectric(1.0) on the NonbondedForce to turn off the reaction field approximation, which does not produce correct results when combined with GBSA.

  • Constructor Details

    • GBSAOBCForce

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

    • addParticle

      public int addParticle(double charge, double radius, double scale)
      Add a particle to the force field.
      Parameters:
      charge - The charge of the particle, measured in units of the proton charge.
      radius - The GBSA radius of the particle, measured in nm.
      scale - The OBC scaling parameter for the particle.
      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 (in nm) being used for nonbonded interactions.
      Returns:
      The cutoff distance, measured in nm.
    • getNonbondedMethod

      public int getNonbondedMethod()
      Get the method used for handling long range nonbonded interactions.
      Returns:
      The nonbonded method.
    • getNumParticles

      public int getNumParticles()
      Get the number of particles in the force field.
      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 scale)
      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 GBSA radius of the particle, measured in nm (output).
      scale - The OBC scaling parameter for the particle (output).
    • getParticleParameters

      public void getParticleParameters(int index, DoubleBuffer charge, DoubleBuffer radius, DoubleBuffer scale)
      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 GBSA radius of the particle, measured in nm (output).
      scale - The OBC scaling parameter for the particle (output).
    • getSoluteDielectric

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

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

      public double getSurfaceAreaEnergy()
      Get the energy scale for the surface area term, measured in kJ/mol/nmˆ2.
      Returns:
      The surface area energy scale.
    • setCutoffDistance

      public void setCutoffDistance(double distance)
      Set the cutoff distance (in nm) being used for nonbonded interactions.
      Parameters:
      distance - The cutoff distance, measured in nm.
    • setNonbondedMethod

      public void setNonbondedMethod(int method)
      Set the method used for handling long range nonbonded interactions.
      Parameters:
      method - The nonbonded method.
    • setParticleParameters

      public void setParticleParameters(int index, double charge, double radius, double scale)
      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 GBSA radius of the particle, measured in nm.
      scale - The OBC scaling parameter for the particle.
    • setSoluteDielectric

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

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

      public void setSurfaceAreaEnergy(double energy)
      Set the energy scale for the surface area term, measured in kJ/mol/nmˆ2.
      Parameters:
      energy - The surface area energy scale.
    • updateParametersInContext

      public void updateParametersInContext(Context context)
      Update the particle 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.