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.
         "translate"    node translates freely - no rotation.
         "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
       
publicIVM.PublicIVM(builtins.object)
IVM

 
class IVM(publicIVM.PublicIVM)
    IVM(sim=0)
 
Internal Variable Module class
members:
  reuse: reuse coordinates and topology from previous call to 
          run(), init(), or autoTorsion()
         this can be unset explicitly, or whenever a topology
         element is modified.
 
 
Method resolution order:
IVM
publicIVM.PublicIVM
builtins.object

Methods defined here:
__init__(s, sim=0)
Constructor. The optional argument is a simulation.Simulation. By
default, the current simulation is used.
addConfigAction(s, action)
Add a command to be run before topology configuration via
protocol.torsionTopology or protocol.cartesianTopology.
autoTorsion(s)
set up hinges for torsion angle dynamics/minimization
baseAtoms(s)
breakAllBondsIn(s, sel)
break all bonds between atoms in selection
breakAllBondsTo(s, sel)
break all bonds to bonded atoms in selection
breakBond(s, *args)
Argument is either two selections specifying one atom each or
a single selection specifying two atoms corresponding to a 
bond to break (topologically). 
The BOND potential is not altered.
calcEnergy(s)
calculate energy in internal coordinates, and return the total
energy. This is a safe wrapper around the underlying calcEnergy
constrainBond(s, el)
specify selection consisting of two atoms corresponding to a
bond to constrain
finalTime(s)
fix(s, sel)
fixedAtoms(s)
Return an atomSel.AtomSel containing all atoms fixed in space.
group(s, sel)
groups(s)
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.]
help(s)
hinge(s, hingeType, sel=None, *sels)
add a hinge definition of type hingeType for atoms specified by
sel, an atomSel.AtomSel, string or list of indices. Some
hinge types require additional selections specified as addtional
arguments.
idAtom(s, i)
info(s, what='default')
return a string with info on the IVM settings.
An optional argument specifies information to return:
   default
   integrate    - integration specific info
   minimization - minimization specific info
   topology     - information about the topology settings
init(s)
iters(s)
Actual number of steps taken during dynamics or minimization. It is
reset and updated by run().
minStepsize(s)
numSteps(s)
outputStepInfo(s, step, stepsize, type, totTime)
printInterval(s)
reset(s)
resetCMInterval(s)
resetReuse(s)
reuse(s)
run(s)
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.
setBaseAtoms(s, l)
setBondExclude(s, v)
setConstraintList(s, v)
setFinalTime(s, v)
setGroupList(s, v)
setHingeList(s, v)
setMinStepsize(s, v)
setNumSteps(s, v)
setPrintInterval(s, v)
setResetCMInterval(s, v)
setStepType(s, v)
setStepsize(s, v)
setStepsizeThreshold(s, v)
setTrajectory(s, v)
stepsize(s)
stepsizeThreshold(s)
time(s)
Actual time elapsed during dynamics. It is reset and updated by
run().
trajectory(s)

Methods inherited from publicIVM.PublicIVM:
Ekinetic(self, *args, **kwargs) -> 'float_type'
Epotential(self, *args, **kwargs) -> 'float_type'
Etotal(self, *args, **kwargs) -> 'float_type'
__repr__ = _swig_repr(self)
adjustStepsize(self, *args, **kwargs) -> 'bool'
bathTemp(self, *args, **kwargs) -> 'float_type'
bondExclude(self, *args, **kwargs) -> 'CDSList< BondIDPair >'
cTolerance(self, *args, **kwargs) -> 'float_type'
calcTemperature(self, *args, **kwargs) -> 'void'
clearRecalc(self, *args, **kwargs) -> 'void'
constrainLengths(self, *args, **kwargs) -> 'bool'
constraintList(self, *args, **kwargs) -> 'CDSList< BondIDPair >'
currentTemp(self, *args, **kwargs) -> 'float_type'
dEpred(self, *args, **kwargs) -> 'float_type'
dim(self, *args, **kwargs) -> 'int'
dof(self, *args, **kwargs) -> 'int'
eCount(self, *args, **kwargs) -> 'int'
eCountReset(self, *args, **kwargs) -> 'void'
eTolerance(self, *args, **kwargs) -> 'float_type'
forceRecalc(self, *args, **kwargs) -> 'void'
frictionCoeff(self, *args, **kwargs) -> 'float_type'
gTolerance(self, *args, **kwargs) -> 'float_type'
gradMagnitude(self, *args, **kwargs) -> 'float_type'
groupList(self, *args, **kwargs) -> 'CDSList< CDSList< int > >'
groupTorsion(self, *args, **kwargs) -> 'void'
hingeList(self, *args, **kwargs) -> 'CDSList< HingeSpec >'
initDynamics(self, *args, **kwargs) -> 'void'
kBoltzmann(self, *args, **kwargs) -> 'float_type'
maxCalls(self, *args, **kwargs) -> 'int'
maxDeltaE(self, *args, **kwargs) -> 'float_type'
maxTSFactor(self, *args, **kwargs) -> 'float_type'
minStepSize(self, *args, **kwargs) -> 'float_type'
minimization(self, *args, **kwargs) -> 'bool'
nodeList(self, *args, **kwargs) -> 'CDSList< PublicNode >'
oldBaseAtoms(self, *args, **kwargs) -> 'CDSList< int >'
pos(self, *args, **kwargs) -> 'CDSVector< double >'
potList(self, *args) -> 'rc_DerivedPot< PotList > &'
printStepDetails(self, *args, **kwargs) -> 'void'
recenterLargeDispl(self, *args, **kwargs) -> 'bool const'
resetCM(self, *args, **kwargs) -> 'void'
responseTime(self, *args, **kwargs) -> 'float_type'
restoreState(self, *args, **kwargs) -> 'void'
saveState(self, *args, **kwargs) -> 'void'
scaleVel(self, *args, **kwargs) -> 'bool'
setAdjustStepsize(self, *args, **kwargs) -> 'void'
setBathTemp(self, *args, **kwargs) -> 'void'
setCTolerance(self, *args, **kwargs) -> 'void'
setConstrainLengths(self, *args, **kwargs) -> 'void'
setDEpred(self, *args, **kwargs) -> 'void'
setETolerance(self, *args, **kwargs) -> 'void'
setFrictionCoeff(self, *args, **kwargs) -> 'void'
setGTolerance(self, *args, **kwargs) -> 'void'
setMaxCalls(self, *args, **kwargs) -> 'void'
setMaxDeltaE(self, *args, **kwargs) -> 'void'
setMaxTSFactor(self, *args, **kwargs) -> 'void'
setMinStepSize(self, *args, **kwargs) -> 'void'
setOldBaseAtoms(self, *args, **kwargs) -> 'void'
setPos(self, *args, **kwargs) -> 'void'
setPotList(self, *args, **kwargs) -> 'void'
setRVecProd(self, *args, **kwargs) -> 'void'
setRVecSize(self, *args, **kwargs) -> 'void'
setRecenterLargeDispl(self, *args, **kwargs) -> 'void'
setResponseTime(self, *args, **kwargs) -> 'void'
setScaleVel(self, *args, **kwargs) -> 'void'
setVecVec3Prod(self, *args, **kwargs) -> 'void'
setVecVec3Size(self, *args, **kwargs) -> 'void'
setVel(self, *args, **kwargs) -> 'void'
setVerbose(self, *args, **kwargs) -> 'void'
simulation(self, *args, **kwargs) -> 'Simulation *'
simulationGroupList(self, *args, **kwargs) -> 'CDSList< CDSList< int > >'
step(self, *args, **kwargs) -> 'int'
stepType(self, *args, **kwargs) -> 'char const *'
vel(self, *args, **kwargs) -> 'CDSVector< double >'
velFromCartesian(self, *args, **kwargs) -> 'void'
verbose(self, *args, **kwargs) -> 'int'

Static methods inherited from publicIVM.PublicIVM:
__swig_destroy__ = delete_PublicIVM(...)

Data descriptors inherited from publicIVM.PublicIVM:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)
thisown

 
The membership flag

Data and other attributes inherited from publicIVM.PublicIVM:
printCMVel = 32
printCoords = 1
printEnergy = 16
printLoopDebug = 4096
printLoopInfo = 8192
printNodeDef = 2048
printNodeForce = 64
printNodeKE = 16384
printNodePos = 128
printNodeTheta = 256
printResetCM = 2
printStepDebug = 512
printStepInfo = 1024
printTemperature = 8
printVelFromCartCost = 4

 
Functions
       
help()
minimizeEscape(ivm)
pyXplorHelp = help()

 
Data
        env = environ({'ARCH': 'Linux_x86_64', 'MWWM': 'allwm'...home/schwitrs/lib/Linux_x86_64::/usr/local/lib'})
minimizeEscapeCnt = 0