Class OpenMMSystem

java.lang.Object
ffx.openmm.System
ffx.potential.openmm.OpenMMSystem
Direct Known Subclasses:
OpenMMDualTopologySystem

public class OpenMMSystem extends System
Create and manage an OpenMM System.

The definition of a System involves four elements:

The particles and constraints are defined directly by the System object, while forces are defined by objects that extend the Force class. After creating a System, call addParticle() once for each particle, addConstraint() for each constraint, and addForce() for each Force.

In addition, particles may be designated as "virtual sites". These are particles whose positions are computed automatically based on the positions of other particles. To define a virtual site, call setVirtualSite(), passing in a VirtualSite object that defines the rules for computing its position.

  • Field Details

    • forceField

      protected ForceField forceField
      The Force Field in use.
    • atoms

      protected Atom[] atoms
      Array of atoms in the system.
    • updateBondedTerms

      protected boolean updateBondedTerms
      This flag indicates bonded force constants and equilibria are updated (e.g. during ManyBody titration).
    • softcoreCreated

      protected boolean softcoreCreated
      Fixed charge softcore vdW force boolean.
    • bondForce

      protected BondForce bondForce
      OpenMM Custom Bond Force
    • angleForce

      protected AngleForce angleForce
      OpenMM Custom Angle Force
    • inPlaneAngleForce

      protected InPlaneAngleForce inPlaneAngleForce
      OpenMM Custom In-Plane Angle Force
    • stretchBendForce

      protected StretchBendForce stretchBendForce
      OpenMM Custom Stretch-Bend Force
    • ureyBradleyForce

      protected UreyBradleyForce ureyBradleyForce
      OpenMM Custom Urey-Bradley Force
    • outOfPlaneBendForce

      protected OutOfPlaneBendForce outOfPlaneBendForce
      OpenMM Custom Out-of-Plane Bend Force
    • piOrbitalTorsionForce

      protected PiOrbitalTorsionForce piOrbitalTorsionForce
      OpenMM Custom Pi-Torsion Force
    • torsionForce

      protected TorsionForce torsionForce
      OpenMM AMOEBA Torsion Force.
    • improperTorsionForce

      protected ImproperTorsionForce improperTorsionForce
      OpenMM Improper Torsion Force.
    • amoebaTorsionTorsionForce

      protected AmoebaTorsionTorsionForce amoebaTorsionTorsionForce
      OpenMM Torsion-Torsion Force.
    • stretchTorsionForce

      protected StretchTorsionForce stretchTorsionForce
      OpenMM Stretch-Torsion Force.
    • angleTorsionForce

      protected AngleTorsionForce angleTorsionForce
      OpenMM Angle-Torsion Force.
    • restrainTorsionsForce

      protected RestrainTorsionsForce restrainTorsionsForce
      OpenMM Restraint-Torsion Force.
    • restrainPositionsForce

      protected RestrainPositionsForce restrainPositionsForce
      OpenMM Restrain-Position Force.
    • restrainGroupsForce

      protected RestrainGroupsForce restrainGroupsForce
      OpenMM Restrain-Groups Force.
    • amoebaVDWForce

      protected AmoebaVdwForce amoebaVDWForce
      OpenMM AMOEBA van der Waals Force.
    • amoebaMultipoleForce

      protected AmoebaMultipoleForce amoebaMultipoleForce
      OpenMM AMOEBA Multipole Force.
    • amoebaGeneralizedKirkwoodForce

      protected AmoebaGeneralizedKirkwoodForce amoebaGeneralizedKirkwoodForce
      OpenMM Generalized Kirkwood Force.
    • amoebaWcaDispersionForce

      protected AmoebaWcaDispersionForce amoebaWcaDispersionForce
      OpenMM AMOEBA WCA Dispersion Force.
    • amoebaGKCavitationForce

      protected AmoebaGKCavitationForce amoebaGKCavitationForce
      OpenMM AMOEBA WCA Cavitation Force.
    • fixedChargeGBForce

      protected FixedChargeGBForce fixedChargeGBForce
      OpenMM Custom GB Force.
    • fixedChargeNonBondedForce

      protected FixedChargeNonbondedForce fixedChargeNonBondedForce
      OpenMM Fixed Charge Non-Bonded Force.
    • fixedChargeAlchemicalForces

      protected FixedChargeAlchemicalForces fixedChargeAlchemicalForces
      Custom forces to handle alchemical transformations for fixed charge systems.
    • andersenThermostat

      protected AndersenThermostat andersenThermostat
      OpenMM thermostat. Currently, an Andersen thermostat is supported.
    • monteCarloBarostat

      protected MonteCarloBarostat monteCarloBarostat
      Barostat to be added if NPT (isothermal-isobaric) dynamics is requested.
    • cmMotionRemover

      protected CMMotionRemover cmMotionRemover
      OpenMM center-of-mass motion remover.
  • Constructor Details

    • OpenMMSystem

      public OpenMMSystem()
      OpenMMSystem constructor.
    • OpenMMSystem

      public OpenMMSystem(OpenMMEnergy openMMEnergy)
      OpenMMSystem constructor.
      Parameters:
      openMMEnergy - ForceFieldEnergyOpenMM instance.
  • Method Details

    • getPotential

      public Potential getPotential()
      Get the Potential in use.
      Returns:
      The Potential.
    • getAtoms

      public Atom[] getAtoms()
      Get the atoms in the system.
      Returns:
      Array of atoms in the system.
    • addForces

      public void addForces()
      Add forces to the system.
    • removeForce

      public void removeForce(Force force)
      Remove a force from the OpenMM System. The OpenMM memory associated with the removed Force object is deleted.
    • addAndersenThermostatForce

      public void addAndersenThermostatForce(double targetTemp)
      Add an Andersen thermostat to the system.
      Parameters:
      targetTemp - Target temperature in Kelvins.
    • addAndersenThermostatForce

      public void addAndersenThermostatForce(double targetTemp, double collisionFreq)
      Add an Andersen thermostat to the system.
      Parameters:
      targetTemp - Target temperature in Kelvins.
      collisionFreq - Collision frequency in 1/psec.
    • addCOMMRemoverForce

      public void addCOMMRemoverForce()
      Adds a force that removes center-of-mass motion.
    • addMonteCarloBarostatForce

      public void addMonteCarloBarostatForce(double targetPressure, double targetTemp, int frequency)
      Add a Monte Carlo Barostat to the system.
      Parameters:
      targetPressure - The target pressure (in atm).
      targetTemp - The target temperature.
      frequency - The frequency to apply the barostat.
    • calculateDegreesOfFreedom

      public int calculateDegreesOfFreedom()
      Calculate the number of degrees of freedom.
      Returns:
      Number of degrees of freedom.
    • getTemperature

      public double getTemperature(double kineticEnergy)
    • getForceField

      public ForceField getForceField()
      Get the ForceField in use.
      Returns:
      ForceField instance.
    • getCrystal

      public Crystal getCrystal()
      Get the Crystal instance.
      Returns:
      the Crystal instance.
    • getNumberOfVariables

      public int getNumberOfVariables()
      Get the number of variables.
      Returns:
      the number of variables.
    • free

      public void free()
      Destroy the system.
    • setUpdateBondedTerms

      public void setUpdateBondedTerms(boolean updateBondedTerms)
    • setDefaultPeriodicBoxVectors

      protected void setDefaultPeriodicBoxVectors()
      Set the default values of the vectors defining the axes of the periodic box (measured in nm).

      Any newly created Context will have its box vectors set to these. They will affect any Force added to the System that uses periodic boundary conditions.

      Triclinic boxes are supported, but the vectors must satisfy certain requirements. In particular, a must point in the x direction, b must point "mostly" in the y direction, and c must point "mostly" in the z direction. See the documentation for details.

    • updateParameters

      public void updateParameters(@Nullable Atom[] atoms)
      Update parameters if the Use flags and/or Lambda value has changed.
      Parameters:
      atoms - Atoms in this list are considered.
    • addAtoms

      protected void addAtoms() throws Exception
      Adds atoms from the molecular assembly to the OpenMM System and reports to the user the number of particles added.
      Throws:
      Exception
    • updateAtomMass

      public boolean updateAtomMass()
      This method sets the mass of inactive atoms to zero.
      Returns:
      Returns true if any inactive atoms were found and set to zero mass.
    • hasAmoebaCavitationForce

      public boolean hasAmoebaCavitationForce()
    • addUpBondConstraints

      protected void addUpBondConstraints()
      Add a constraint to every bond.
    • addHydrogenConstraints

      protected void addHydrogenConstraints()
      Add a constraint to every bond that includes a hydrogen atom.
    • setUpHydrogenAngleConstraints

      protected void setUpHydrogenAngleConstraints()
      Add a constraint to every angle that includes two hydrogen atoms.
    • isHydrogenAngle

      protected boolean isHydrogenAngle(Angle angle)
      Check to see if an angle is a hydrogen angle. This method only returns true for hydrogen angles that are less than 160 degrees.
      Parameters:
      angle - Angle to check.
      Returns:
      boolean indicating whether an angle is a hydrogen angle that is less than 160 degrees.