ivm |
index |
The Internal Variable Module
An IVM object manipulates atomic coordinates, using molecular
dynamics and gradient minimization techniques. IVM objects efficiently
perform calculations in general coordinates. Particularly useful are
calculations in torsion angle space, and also rigid-body
manipulations.
A copy of the IVM paper with some corrections is available online.
The proper citation is
C.D. Schwieters and G.M. Clore, Internal coordinates for
molecular dynamics and minimization in structure determination and
refinement; J. Magn. Reson. 152, 288-302 (2001).
Constructor:
IVM(simulation) - the optional simulation argument specifies the
associated simulation.Simulation and defaults to
currentSimulation if omitted.
Members:
run() - perform dynamics or minimization. Stop when numSteps have been taken,
or finalTime has been reached (whichever is first). Returns an
integer which is 1 on successfull completion of dynamics or
minimization and negative on failure.
iters() - the number of dynamics or minimization steps actually taken
time() - the actual total time elapsed during dynamics.
reset() - clear the following lists:
bondExclude, groupList, hingeList, constraintList, baseAtoms
fix(selList) - add fixed atoms specified by a selection of a
list of indices or atoms. Note that this command modifies
groupList.
group(selList) - add atoms to groupList. Atoms may specified by a
selection of a list of indices or atoms.
hinge(type,specification) - add a hinge specification. The specification
should be a sequence, the first element of which is a
string hinge type and the second a atom selection.
Any optional arguments then follow.
Possible hinge types:
"full" node translates and rotates freely.
"rotate" node rotates freely - no translation.
"rotate2" node rotates in two dimensions. This option has
not been fully tested.
"translate" node translates freely - no rotation.
"translate1" No rotation. Node only translates in the direction
given by atom2.pos() - atom1.pos(), where these
atoms are specified by two additional selection
arguments.
"torsion" allow rotation of child hinge about bond
connecting it to its parent.
"bend" allow only bending motion.
When using this hinge type, sel should specify
only a single atom, and two additonal
selection argument are required. The second selection
should specify the parent atom, and the third
should define the bend angle.
"bendtorsion" in addition to torsion motion, allow
bending. Specifications for the extra atom selection
argument are the same as for bend hinges.
breakBond(a,b) - a sequence of two atom selections. This adds an
entry to the bondExclude list.
breakAllBondsTo(sel) - break all bonds to atoms in the the given selection.
breakAllBondsIn(sel) - break all bonds between the atoms in the
given selection.
constrainBond(sel) - sel should select two bound atoms. The corresponding
bond length will be constrained if
setConstrainLengths(True) is also called.
autoTorsion() - set up hinges for torsion angle dynamics/minimization
groups() - Return a list of atomSel.AtomSel objects corresponding
to atoms grouped in rigid bodies in this IVM object.
[This method converts the indices from groupList()
into AtomSel objects, discarding negative indices
which are used for IVM internal purposes.]
velFromCartesian() - set internal variable velocities to be as close
as possible to the current Cartesian counterparts.
info(what) - return a string with info on the IVM settings.
Optional argument what specifies information to return:
default - info common to integration and minimization
all - everything
integrate - integration specific info
minimization - minimization specific info
topology - information about the topology settings
Accessors:
following are functions which return values of the quantities.
To set the values, use functions of the form setQuantity( value ).
Note that not all quantities are accessible.
stepType() - integration or minimization method. Possible values:
minimization: Powell, MinimizeCG, ConMin, Simplex
integration: PC6, RungeKutta, Gear, Milne, Verlet
dof() - return the number of degrees of freedom
dim() - return the dimension: number of internal coordinates
(this is usually larger than dof() because there are
additional constraint relationships).
pos() - values of the internal coordinates (NOT atom
positions). Note that setPos() also modifies ic
velocities to be consistent with any constraints.
vel() - values of the internal coordinate velocities. Note
that velocities actually set maybe modified so they
areconsistent with any constraints. If stepType is
set to a minimization type, velocities are zero'd
when the run() method is called.
numSteps() - number of time steps to take. 0 for disabled
finalTime() - end time for integration. 0 for disabled
stepsize() - the current step size.
printInterval() - how many steps between energy printouts. Set to 0
for no output.
resetCMInterval() - how many steps between zeroing linear/ang momentum.
autoTorsion() - set up groupList and hingeList for
torsion-angle manipulations. Existing groupList and
and hingeList specifications are preserved.
bondExclude() - list of pairs of inidices:
bonds to ignore when considering molecular topology
constraintList() - list of pairs of inidices:
bonds to constrain (if useLengthConstraints)
groupList() - list of list of atom indices-
grouped atoms fixed w/ respect to each other
simulationGroupList() - the same as groupList() except indices corresponding
to atoms not in the simulation() are omitted.
hingeList() - list of atom indices:
atom(s) to use as base atoms in the tree topology
baseAtoms() - list of atom indices:
atom(s) which should be used as base atoms.
oldBaseAtoms() - list of atom indices:
atom(s) which were base atoms previously
nodeList() - return list of nodes, each of type PublicNode,
described below.
eCount() - number of energy evaluations since last eReset()
bathTemp() - value bath temp is set to
currentTemp() - temperature calc'd from current velocities
Etotal() - current total energy (pot + kin)
Epotential() - current potential energy
Ekinetic() - current kinetic energy
eTolerance() - energy tolerance - for minimization and integration
gTolerance() - gradient tolerance - for minimization
cTolerance() - length constraint tolerance
dEpred() - delta energy prediction
responseTime() - response time for timestep adjustment
frictionCoeff() - friction scaling coefficient
(should be zero if velocity scaling is used)
kBoltzmann() - Boltzmann's constant
maxCalls() - maximum number of energy calls in minimization
constrainLengths() - boolean: whether to use length constraints
adjustStepsize() - boolean: whether to adjust integration stepsize
to maintain constant energy.
scaleVel() - boolean: whether to scale velocities for
equilibration with a temperature bath.
maxTSFactor() - max factor for timestep increase (1.025)
maxDeltaE() - maximum energy change before halving timestep
minStepsize() - stepsize below which an exception is thrown by
step()
stepsizeThreshold() - When the minStepsize limit is reached in dynamics,
gradient minimization is performed if the stepsize
is smaller than this value (default: 0).
trajectory() - the trajectory method. This object must have a
method named write(stepNum,time).
recenterLargeDispl() - recenter all atoms at the origin if any atom
has an espcially large displacement (False).
verbose() - return the verbose bitfield: Possible values:
printCoords
printResetCM
printVelFromCartCost
printTemperature
printEnergy
printCMVel
printNodeForce
printNodePos
printNodeTheta
printStepDebug
printStepInfo
printNodeDef
printLoopDebug
printLoopInfo
an example of using verbose:
min.setVerbose( min.printLoopInfo | min.printStepInfo )
The following methods are primarily for internal use:
idAtom(index) - given an atom index, return descriptive string
init() - initialize IVM internals
initDynamics(reuseTopology) -
For dynamics, resynchronize internal positions
and velocities from Cartesian values. For
minimization, velocities are zeroed.
Also. reinitialize the integrator.
If reuseTopology is false, regenerate topology.
calcEnergy() - energy and derivatives calc'd and terms placed in energyTerms.
The toal energy (pot+kinetic) is returned.
Note that this routine calls ivm.init() so that ic tables are
reset if any topology changes have been made. Etotal,
Epotential and Ekinetic are updated using this call.
eCountReset() - reset eCount, incremented during each energy calc
groupTorsion() - group atoms for torsion-angle dynamics. Obeys pre-existing
groupList, hingeList, bondExclude and constraintList.
resetCM() - zero center of mass linear and angular momentum
step(stepsize) - take dynamics or minimization step
tuple (finished,stepsize) returned.
Finished is true if stop condition is met.
stepsize specifies new stepsize, if variable.
class PublicNode, a list of which is returned by the nodeList() member of IVM.
Methods:
dim() - return the dimension of the hinge. This includes the number
of degrees of freedom plus the number of constraints.
startIndex() - the start index of coordinates for this hinge in the
IVM's pos() and vel() arrays of internal coordinates.
parentAtom() - index of atom in the parent node to which this hinge/node is
attached.
type() - type of hinge.
atoms() - list of indices of atoms contained in this node.
atoms()[0] is the index of the hinge atom.
Classes | ||||||||||
|
Functions | ||
|
Data | ||
env = environ({'ARCH': 'Linux_x86_64', 'MWWM': 'allwm'...ib/Linux_x86_64::/usr/local/lib:/usr/local/lib'}) minimizeEscapeCnt = 0 |