OpenMM
 All Classes Namespaces Functions Variables Pages
CMAPTorsionForce Class Reference

This class implements an interaction between pairs of dihedral angles. More...

+ Inheritance diagram for CMAPTorsionForce:

Public Member Functions

def getNumMaps
 getNumMaps(CMAPTorsionForce self) -> int More...
 
def getNumTorsions
 getNumTorsions(CMAPTorsionForce self) -> int More...
 
def addMap
 addMap(CMAPTorsionForce self, int size, vectord energy) -> int More...
 
def getMapParameters
 getMapParameters(CMAPTorsionForce self, int index) More...
 
def setMapParameters
 setMapParameters(CMAPTorsionForce self, int index, int size, vectord energy) More...
 
def addTorsion
 addTorsion(CMAPTorsionForce self, int map, int a1, int a2, int a3, int a4, int b1, int b2, int b3, int b4) -> int More...
 
def getTorsionParameters
 getTorsionParameters(CMAPTorsionForce self, int index) More...
 
def setTorsionParameters
 setTorsionParameters(CMAPTorsionForce self, int index, int map, int a1, int a2, int a3, int a4, int b1, int b2, int b3, int b4) More...
 
def __init__
 init(OpenMM::CMAPTorsionForce self) -> CMAPTorsionForce init(OpenMM::CMAPTorsionForce self, CMAPTorsionForce other) -> CMAPTorsionForce More...
 
def __del__
 del(OpenMM::CMAPTorsionForce self) More...
 
- Public Member Functions inherited from Force
def __init__
 
def __del__
 del(OpenMM::Force self) More...
 
def getForceGroup
 getForceGroup(Force self) -> int More...
 
def setForceGroup
 setForceGroup(Force self, int group) More...
 
def __copy__
 
def __deepcopy__
 

Public Attributes

 this
 

Detailed Description

This class implements an interaction between pairs of dihedral angles.

The interaction energy is defined by an "energy correction map" (CMAP), which is simply a set of tabulated energy values on a regular grid of (phi, psi) angles. Natural cubic spline interpolation is used to compute forces and energies at arbitrary values of the two angles.

To use this class, first create one or more energy correction maps by calling addMap(). For each one, you provide an array of energies at uniformly spaced values of the two angles. Next, add interactions by calling addTorsion(). For each one, you specify the sequence of particles used to calculate each of the two dihedral angles, and the index of the map used to calculate their interaction energy.

Constructor & Destructor Documentation

def __init__ (   self,
  args 
)

init(OpenMM::CMAPTorsionForce self) -> CMAPTorsionForce init(OpenMM::CMAPTorsionForce self, CMAPTorsionForce other) -> CMAPTorsionForce

Create a CMAPTorsionForce.

References simtk.openmm.openmm.stripUnits().

def __del__ (   self)

del(OpenMM::CMAPTorsionForce self)

References simtk.openmm.openmm.stripUnits().

Member Function Documentation

def addMap (   self,
  args 
)

addMap(CMAPTorsionForce self, int size, vectord energy) -> int

Create a new map that can be used for torsion pairs.

Parameters
sizethe size of the map along each dimension
energythe energy values for the map. This must be of length size*size. The element energy[i+size*j] contains the energy when the first torsion angle equals i*2*PI/size and the second torsion angle equals j*2*PI/size.

References simtk.openmm.openmm.stripUnits().

def addTorsion (   self,
  args 
)

addTorsion(CMAPTorsionForce self, int map, int a1, int a2, int a3, int a4, int b1, int b2, int b3, int b4) -> int

Add a CMAP torsion term to the force field.

Parameters
mapthe index of the map to use for this term
a1the index of the first particle forming the first torsion
a2the index of the second particle forming the first torsion
a3the index of the third particle forming the first torsion
a4the index of the fourth particle forming the first torsion
b1the index of the first particle forming the second torsion
b2the index of the second particle forming the second torsion
b3the index of the third particle forming the second torsion
b4the index of the fourth particle forming the second torsion

References simtk.openmm.openmm.stripUnits().

def getMapParameters (   self,
  args 
)

getMapParameters(CMAPTorsionForce self, int index)

Get the energy values of a map.

Parameters
indexthe index of the map for which to get energy values
sizethe size of the map along each dimension
energythe energy values for the map. This must be of length size*size. The element energy[i+size*j] contains the energy when the first torsion angle equals i*2*PI/size and the second torsion angle equals j*2*PI/size.

References simtk.openmm.openmm.stripUnits().

def getNumMaps (   self)

getNumMaps(CMAPTorsionForce self) -> int

Get the number of maps that have been defined.

References simtk.openmm.openmm.stripUnits().

def getNumTorsions (   self)

getNumTorsions(CMAPTorsionForce self) -> int

Get the number of CMAP torsion terms in the potential function

References simtk.openmm.openmm.stripUnits().

def getTorsionParameters (   self,
  args 
)

getTorsionParameters(CMAPTorsionForce self, int index)

Get the force field parameters for a CMAP torsion term.

Parameters
indexthe index of the torsion for which to get parameters
mapthe index of the map to use for this term
a1the index of the first particle forming the first torsion
a2the index of the second particle forming the first torsion
a3the index of the third particle forming the first torsion
a4the index of the fourth particle forming the first torsion
b1the index of the first particle forming the second torsion
b2the index of the second particle forming the second torsion
b3the index of the third particle forming the second torsion
b4the index of the fourth particle forming the second torsion

References simtk.openmm.openmm.stripUnits().

def setMapParameters (   self,
  args 
)

setMapParameters(CMAPTorsionForce self, int index, int size, vectord energy)

Set the energy values of a map.

Parameters
indexthe index of the map for which to set energy values
sizethe size of the map along each dimension
energythe energy values for the map. This must be of length size*size. The element energy[i+size*j] contains the energy when the first torsion angle equals i*2*PI/size and the second torsion angle equals j*2*PI/size.

References simtk.openmm.openmm.stripUnits().

def setTorsionParameters (   self,
  args 
)

setTorsionParameters(CMAPTorsionForce self, int index, int map, int a1, int a2, int a3, int a4, int b1, int b2, int b3, int b4)

Set the force field parameters for a CMAP torsion term.

Parameters
indexthe index of the torsion for which to set parameters
mapthe index of the map to use for this term
a1the index of the first particle forming the first torsion
a2the index of the second particle forming the first torsion
a3the index of the third particle forming the first torsion
a4the index of the fourth particle forming the first torsion
b1the index of the first particle forming the second torsion
b2the index of the second particle forming the second torsion
b3the index of the third particle forming the second torsion
b4the index of the fourth particle forming the second torsion

References simtk.openmm.openmm.stripUnits().

Member Data Documentation

this

Referenced by System.__init__().


The documentation for this class was generated from the following file: