Package ffx.openmm

Class NoseHooverIntegrator

java.lang.Object
ffx.openmm.Integrator
ffx.openmm.NoseHooverIntegrator
Direct Known Subclasses:
DrudeNoseHooverIntegrator

public class NoseHooverIntegrator extends Integrator
This is an Integrator which simulates a System using one or more Nose Hoover chain thermostats, using the "middle" leapfrog propagation algorithm described in J. Phys. Chem. A 2019, 123, 6056-6079.
  • Constructor Details

    • NoseHooverIntegrator

      public NoseHooverIntegrator(com.sun.jna.ptr.PointerByReference pointer)
      Create a NoseHooverIntegrator with a single thermostat.
      Parameters:
      pointer - A pointer to the native OpenMM NoseHooverIntegrator object.
    • NoseHooverIntegrator

      public NoseHooverIntegrator(double stepSize)
      Create a NoseHooverIntegrator with a single thermostat.
      Parameters:
      stepSize - The step size with which to integrate the system (in ps).
    • NoseHooverIntegrator

      public NoseHooverIntegrator(double stepSize, double temperature, double collisionFrequency, int numMTS, int numYoshidaSuzuki, int numNoseHoover)
      Create a NoseHooverIntegrator with detailed thermostat parameters.
      Parameters:
      stepSize - The step size with which to integrate the system (in ps).
      temperature - The temperature of the heat bath (in Kelvin).
      collisionFrequency - The collision frequency (in 1/ps).
      numMTS - The number of multiple time step levels.
      numYoshidaSuzuki - The number of Yoshida-Suzuki steps.
      numNoseHoover - The number of Nosé-Hoover chain thermostats.
  • Method Details

    • addSubsystemThermostat

      public int addSubsystemThermostat(com.sun.jna.ptr.PointerByReference particles, com.sun.jna.ptr.PointerByReference chainWeights, double temperature, double collisionFrequency, double relativeTemperature, double relativeCollisionFrequency, int numMTS, int numYoshidaSuzuki, int numNoseHoover)
      Add a subsystem thermostat to the integrator.
      Parameters:
      particles - The particles to be controlled by this thermostat.
      chainWeights - The weights for the thermostat chain.
      temperature - The temperature for this thermostat (in Kelvin).
      collisionFrequency - The collision frequency for this thermostat (in 1/ps).
      relativeTemperature - The relative temperature scaling factor.
      relativeCollisionFrequency - The relative collision frequency scaling factor.
      numMTS - The number of multiple time step levels.
      numYoshidaSuzuki - The number of Yoshida-Suzuki steps.
      numNoseHoover - The number of Nosé-Hoover chain thermostats.
      Returns:
      The index of the thermostat that was added.
    • addThermostat

      public int addThermostat(double temperature, double collisionFrequency, int numMTS, int numYoshidaSuzuki, int numNoseHoover)
      Add a thermostat to the integrator.
      Parameters:
      temperature - The temperature for this thermostat (in Kelvin).
      collisionFrequency - The collision frequency for this thermostat (in 1/ps).
      numMTS - The number of multiple time step levels.
      numYoshidaSuzuki - The number of Yoshida-Suzuki steps.
      numNoseHoover - The number of Nosé-Hoover chain thermostats.
      Returns:
      The index of the thermostat that was added.
    • computeHeatBathEnergy

      public double computeHeatBathEnergy()
      Compute the total energy of all heat baths.
      Returns:
      The total heat bath energy.
    • destroy

      public void destroy()
      Destroy the integrator.
      Specified by:
      destroy in class Integrator
    • getCollisionFrequency

      public double getCollisionFrequency(int thermostat)
      Get the collision frequency for a thermostat.
      Parameters:
      thermostat - The index of the thermostat.
      Returns:
      The collision frequency (in 1/ps).
    • getMaximumPairDistance

      public double getMaximumPairDistance()
      Get the maximum pair distance for neighbor list updates.
      Returns:
      The maximum pair distance (in nm).
    • getNumThermostats

      public int getNumThermostats()
      Get the number of thermostats.
      Returns:
      The number of thermostats.
    • getRelativeCollisionFrequency

      public double getRelativeCollisionFrequency(int thermostat)
      Get the relative collision frequency for a thermostat.
      Parameters:
      thermostat - The index of the thermostat.
      Returns:
      The relative collision frequency scaling factor.
    • getRelativeTemperature

      public double getRelativeTemperature(int thermostat)
      Get the relative temperature for a thermostat.
      Parameters:
      thermostat - The index of the thermostat.
      Returns:
      The relative temperature scaling factor.
    • getTemperature

      public double getTemperature(int thermostat)
      Get the temperature for a thermostat.
      Parameters:
      thermostat - The index of the thermostat.
      Returns:
      The temperature (in Kelvin).
    • getThermostat

      public com.sun.jna.ptr.PointerByReference getThermostat(int thermostat)
      Get a reference to a thermostat.
      Parameters:
      thermostat - The index of the thermostat.
      Returns:
      A reference to the thermostat object.
    • hasSubsystemThermostats

      public int hasSubsystemThermostats()
      Check if the integrator has subsystem thermostats.
      Returns:
      1 if subsystem thermostats are present, 0 otherwise.
    • setCollisionFrequency

      public void setCollisionFrequency(double collisionFrequency, int thermostat)
      Set the collision frequency for a thermostat.
      Parameters:
      collisionFrequency - The collision frequency (in 1/ps).
      thermostat - The index of the thermostat.
    • setMaximumPairDistance

      public void setMaximumPairDistance(double distance)
      Set the maximum pair distance for neighbor list updates.
      Parameters:
      distance - The maximum pair distance (in nm).
    • setRelativeCollisionFrequency

      public void setRelativeCollisionFrequency(double relativeCollisionFrequency, int thermostat)
      Set the relative collision frequency for a thermostat.
      Parameters:
      relativeCollisionFrequency - The relative collision frequency scaling factor.
      thermostat - The index of the thermostat.
    • setRelativeTemperature

      public void setRelativeTemperature(double relativeTemperature, int thermostat)
      Set the relative temperature for a thermostat.
      Parameters:
      relativeTemperature - The relative temperature scaling factor.
      thermostat - The index of the thermostat.
    • setTemperature

      public void setTemperature(double temperature, int thermostat)
      Set the temperature for a thermostat.
      Parameters:
      temperature - The temperature (in Kelvin).
      thermostat - The index of the thermostat.
    • step

      public void step(int steps)
      Advance a simulation through time by taking a series of time steps.
      Overrides:
      step in class Integrator
      Parameters:
      steps - The number of time steps to take.