OpenMM
AmberPrmtopFile Class Reference

AmberPrmtopFile parses an AMBER prmtop file and constructs a Topology and (optionally) an OpenMM System from it. More...

List of all members.

Public Member Functions

def __init__
 Load a prmtop file.
def createSystem
 Construct an OpenMM System representing the topology described by this prmtop file.

Public Attributes

 topology
 The Topology read from the prmtop file.
 elements

Detailed Description

AmberPrmtopFile parses an AMBER prmtop file and constructs a Topology and (optionally) an OpenMM System from it.


Constructor & Destructor Documentation

def __init__ (   self,
  file 
)

Load a prmtop file.


Member Function Documentation

def createSystem (   self,
  nonbondedMethod = ff.NoCutoff,
  nonbondedCutoff = 1.0*u.nanometer,
  constraints = None,
  rigidWater = True,
  implicitSolvent = None,
  implicitSolventSaltConc = 0.0*(u.moles/u.liter),
  implicitSolventKappa = None,
  temperature = 298.15*u.kelvin,
  soluteDielectric = 1.0,
  solventDielectric = 78.5,
  removeCMMotion = True,
  hydrogenMass = None,
  ewaldErrorTolerance = 0.0005,
  switchDistance = 0.0*u.nanometer 
)

Construct an OpenMM System representing the topology described by this prmtop file.

Parameters:
nonbondedMethod(object=NoCutoff) The method to use for nonbonded interactions. Allowed values are NoCutoff, CutoffNonPeriodic, CutoffPeriodic, Ewald, or PME.
nonbondedCutoff(distance=1*nanometer) The cutoff distance to use for nonbonded interactions
constraints(object=None) Specifies which bonds angles should be implemented with constraints. Allowed values are None, HBonds, AllBonds, or HAngles.
rigidWater(boolean=True) If true, water molecules will be fully rigid regardless of the value passed for the constraints argument
implicitSolvent(object=None) If not None, the implicit solvent model to use. Allowed values are HCT, OBC1, OBC2, GBn, or GBn2.
implicitSolventSaltConc(float=0.0*unit.moles/unit.liter) The salt concentration for GB calculations (modelled as a debye screening parameter). It is converted to the debye length (kappa) using the provided temperature and solventDielectric
temperature(float=300*kelvin) Temperature of the system. Only used to compute the Debye length from implicitSolventSoltConc
implicitSolventKappa(float units of 1/length) If this value is set, implicitSolventSaltConc will be ignored.
soluteDielectric(float=1.0) The solute dielectric constant to use in the implicit solvent model.
solventDielectric(float=78.5) The solvent dielectric constant to use in the implicit solvent model.
removeCMMotion(boolean=True) If true, a CMMotionRemover will be added to the System
hydrogenMass(mass=None) The mass to use for hydrogen atoms bound to heavy atoms. Any mass added to a hydrogen is subtracted from the heavy atom to keep their total mass the same.
ewaldErrorTolerance(float=0.0005) The error tolerance to use if nonbondedMethod is Ewald or PME.
switchDistance(float=0*nanometers) The distance at which the potential energy switching function is turned on for Lennard-Jones interactions. If the switchDistance is 0 or evaluates to boolean False, no switching function will be used. Values greater than nonbondedCutoff or less than 0 raise ValueError
Returns:
(System) the newly created System

Member Data Documentation

The Topology read from the prmtop file.


The documentation for this class was generated from the following file:
 All Classes Functions Variables