protocol
index


 
High-level initialization routines.
 
These routines provide high level interfaces to initialize potential terms
and to set up minimization and dynamics calculations.
 
The following functions are imported from the regularize module:
fixupCovalentGeom, covalentMinimize, addUnknownAtoms.
 
Global variable pdbLocation: location from which PDB entries are loaded.

 
Classes
       
builtins.object
InitCoordsResult
builtins.tuple(builtins.object)
TopologyEntry

 
class InitCoordsResult(builtins.object)
    A class to simplify returning B factors, occupancies, remarks, and 
unrecognized ATOM records read members, bfactors and occupancies contain 
the respective properties from the PDB record for each atom references by 
its index. ATOM (or optionally HETATM) records which can't be matched to
the loaded PSF are contained in the unknownEntries member.
 
  Methods defined here:
__init__(self)
Initialize self.  See help(type(self)) for accurate signature.

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)

 
class TopologyEntry(builtins.tuple)
    TopologyEntry(filename, linkages, firstStatements, lastStatements, options)
 
TopologyEntry(filename, linkages, firstStatements, lastStatements, options)
 
 
Method resolution order:
TopologyEntry
builtins.tuple
builtins.object

Methods defined here:
__getnewargs__(self)
Return self as a plain tuple.  Used by copy and pickle.
__repr__(self)
Return a nicely formatted representation string
_asdict(self)
Return a new dict which maps field names to their values.
_replace(self, /, **kwds)
Return a new TopologyEntry object replacing specified fields with new values

Class methods defined here:
_make(iterable) from builtins.type
Make a new TopologyEntry object from a sequence or iterable

Static methods defined here:
__new__(_cls, filename, linkages, firstStatements, lastStatements, options)
Create new instance of TopologyEntry(filename, linkages, firstStatements, lastStatements, options)

Data descriptors defined here:
filename

 
Alias for field number 0
linkages

 
Alias for field number 1
firstStatements

 
Alias for field number 2
lastStatements

 
Alias for field number 3
options

 
Alias for field number 4

Data and other attributes defined here:
__match_args__ = ('filename', 'linkages', 'firstStatements', 'lastStatements', 'options')
_field_defaults = {}
_fields = ('filename', 'linkages', 'firstStatements', 'lastStatements', 'options')

Methods inherited from builtins.tuple:
__add__(self, value, /)
Return self+value.
__contains__(self, key, /)
Return key in self.
__eq__(self, value, /)
Return self==value.
__ge__(self, value, /)
Return self>=value.
__getattribute__(self, name, /)
Return getattr(self, name).
__getitem__(self, key, /)
Return self[key].
__gt__(self, value, /)
Return self>value.
__hash__(self, /)
Return hash(self).
__iter__(self, /)
Implement iter(self).
__le__(self, value, /)
Return self<=value.
__len__(self, /)
Return len(self).
__lt__(self, value, /)
Return self<value.
__mul__(self, value, /)
Return self*value.
__ne__(self, value, /)
Return self!=value.
__rmul__(self, value, /)
Return value*self.
count(self, value, /)
Return number of occurrences of value.
index(self, value, start=0, stop=9223372036854775807, /)
Return first index of value.
 
Raises ValueError if the value is not present.

Class methods inherited from builtins.tuple:
__class_getitem__(...) from builtins.type
See PEP 585

 
Functions
       
addDisulfideBond(sel1, sel2)
deprecated. Calls psfGen.addDisulfideBond
addMassSetup(massFunc, name)
register a mass setup function to be called by the massSetup
helper function. The second argument should identify the module in which
the mass setup function is defined.
addPseudoResName(resname)
Add a residue name to the set which correspond to pseudo atoms.
 
This function only updates the list of residue names which correspond
to pseudoatoms. The function updatePseudoAtoms must be called after atoms
are created to actually update the PSEUdo selection.
addTopologySetup(topologyFunc, name)
register a topology setup function to be called by the
torsionTopology and cartesianTopology helper functions. The second
argument should identify the module in which the topology setup
function is defined.
cartesianTopology(ivm, sel='not PSEUDO', oTensors=[])
configure the ivm.IVM tolopogy for Cartesian dynamics/minimization -
for the specified selection.
 
This consists of breaking topological bonds, and specifying that all atoms
are tree ``bases.''
 
This function should be called after any custom changes are made to the
ivm's topology setup, but before torsionTopology().
 
oTensors is a list of varTensor.VarTensor objects used for alignment
tensors.
checkPseudoTopology(ivm)
Check that all pseudo atoms are explicitly configured in the specified 
ivm.IVM object. Otherwise, print a warning.
commandNoEcho(xSim, cmd, verbose=None)
run the given XPLOR command in the specified
xplorSimulation.XplorSimulation.
downloadPDB(entry, content='pdb')
download a compressed PDB file from the wwPDB database
return the PDB record as a string. If the download fails, the function
will sleep for three seconds and then retry. The total number of connection
attempts is set via the module scope variable downloadPDB_connectAttempts
whose default value is 3.
 
The content argument can also be set to 'cifHeader', in which case only their
mmCif header of the specified PDB entry is returned.
genExtendedStructure(pdbFilename=0, selection='not pseudo', verbose=0, maxFixupIters=500, sel=0)
Generate an arbitrary extended conformation with ideal covalent geometry.
 
This assigns X, Y, and Z coordinates to each atom, and then corrects the
covalent geometry with protocol.fixupCovalentGeom(), using maxFixupIters
(an integer) as the maxIters argument to that function (maxFixupIters=0
results in non-ideal covalent geometry).  Note that regardless of the
choice of maxFixupIters, the nonbonded interactions may not be totally
satisfied.
 
The Y and Z coordinates are random (but small enough (within a range of -0.5
to 0.5) to allow bonded atoms to form their bonds) and the X coordinate is
the atom number divided by 10.  This will result in an extended
configuration along the X axis.
 
If a string is provided as pdbFilename, the generated coordinates will be
written to a pdb file named pdbFilename.  (Note that if a pdb file with that
name already exists, it will be used by protocol.fixupCovalentGeom(), 
and the resulting coordinates will update the input file.)
 
The selection argument can be used to specify a subset of atoms to act on. 
This argument replaces the deprecated sel argument.
genTopParFilename(filename, useTopParVar=False, suffix='')
If the given filename exists, it is returned, else the file is
looked for the in the TOPPAR directory, and a full path is returned
if it is found there.
 
If useTopParVar=True, the XPLOR TOPPAR variable is used for the
filename (instead of expanding its value). This results in
shorter paths, but the path is only valid within the XPLOR interpreter.
 
If suffix is specified, the filename will be tried with "." + suffix
appended.
initBiomtCoords(biomt, simulation=None)
Fill in coordinates from a psfGen.Biomt object.
 
Note that occupancy and bfactor values for duplicated coordinates will
not be filled.
initCarb13(filenames=[], scale=0.5)
Initialize the XPLOR carbon chemical shift potential term.
 
arguments are:
   filenames   - either a single filename, or a sequence of filenames of
                 chemical shift assignment tables.
 
In addition to these files C13SHIFTS:rcoil_c13.tbl and
C13SHIFTS:expected_edited.tbl are always added.
initCollapse(sel='all', Rtarget=None, scale=1)
Initialize the XPLOR radius of gyration potential term.
 
Parameters are:
 
   sel - Selection string or atomSel.AtomSel instance specifying atoms to
         include in the calculation of the radius of gyration (Rgyr)
         (default: all atoms).
         
   Rtarget - Target Rgyr.  In practice the target Rgyr should be smaller
             than the expected value to force compaction.  If Rtarget is
             not specified, it is predicted from the number of residues, N,
             by
                Rtarget = 2.2 * N^0.38 - 1   (in Angstroms)
             (i.e., the expected Rgyr minus one).
             
   scale - Scale factor (defaults to 1).  Note that the per-assignment scale
           is always 100, so the energy is actually scaled by 100*scale.
initCoords(files=[], string='', model=-1, verbose=-1, erase=False, useChainID=True, includeHETATM=False, correctSymmetricSidechains=True, fixMethylImpropers=False, fixProtonImpropers=False, strictResNames=False, maxUnreadEntries=(20, 0.2), deleteUnknownAtoms=False, loadOccupanciesBfactors=False, allowedAltLoc=None, biomt=None, selection='all', **kwargs)
Initialize coordinates from one of more pdb file,
  or from a string containing a PDB entry.
 
If model is specified, the specified PDB MODEL record will be read.
 
The selection argument should be a string or atomSel.AtomSel object 
specifying which atom's coordinates to be read. It defaults to the string 
"all" corresponding to all atoms in the current Simulation.
 
If erase is set to True, then atom positions in selection are cleared
before reading the pdb record.
 
If useChainID is True (the default), then the chain ID PDB field is
used as the segment name (if set), If the chain ID field is blank it is
ignored regardless.
 
By default, HETATM PDB records are ignored. If includeHETATM=True, they
will be read and treated as ATOM records.
 
If deleteUnknownAtoms is True, then atoms whose coordinates are still
unknown after the PDB file is read will be deleted.
 
If correctSymmetricSidechains is True (the default), then the sidechains
of phe, tyr, asp, and glu residues are checked to ensure that their
symmetric torsion angles (chi3 for glu, chi2 for the rest) are in
the correct (0..180 degree) range. If they are incorrect, then atomic
positions are swapped using selectTools.correctSymmetricSidechains.
 
If strictResNames is True, then residue names in the PDB must match
those in the PSF, otherwise, the residue names need not match
and backbone atoms for mutants can be read. Note that sidechain atoms (or
base atoms in nucleic acids) may be strangely placed if this option is
False.
 
PDB ATOM records which do not exactly match the expected input
(the atom name is not exactly the same) are attempted to be
matched using the function matchInexactAtomEntry. The argument
maxUnreadEntries specifies the behavior if many unmatachable
records are encountered. If maxUnreadEntries is an integer, an
exception is thrown if more than that number of unmatchable
records are encountered. If maxUnreadEntries is a two-membered
tuple, the exception is not thrown if number of unmatchable
records is less than the first entry or the number of unmatchable
records is less than the second entry (a fraction) multiplied by
the number of atoms. This check is disabled (no exception will be
thrown) if maxUnreadEntries is set to None. The check for unreadable 
entries is also disabled if the selection argument is set to a value
other than "all".
 
If loadOccupanciesBfactors is True, occupancies and bfactors will be
passed to the appropriate xplorSimulation.XplorSimulation. If
specified as True, the atomSelLang facility will support searching
the b and q attritbutes.
 
biomt is an optional psfGen.Biomt object used to generate identical
copies of coordinates in ATOM records.
 
allowedAltLoc specifies the subset of ATOM records (and HETATM records, 
if includeHETATM=True) are read by the contents of the altLoc field. By 
default, all records with a space or the character 'A' in the altLoc field
are read.
 
If fixMethylImpropers is True, selectTools.fixMethylImpropers will be
called after coordinates are loaded.
 
If fixProtonImpropers is True, selectTools.fixProtonImpropers will be
called after coordinates are loaded.
 
swapProtons and assumeSegid are handled by the function
matchInexactAtomEntry. 
 
This function results an InitCoordsResult structure containing
bfactors, occupancies, remarks, and numAtomsRead from the PDB record. 
The numAtomsRead member corresponds to the number read using pdbTool plus 
the entries added with matchInexactCoords.
initDihedrals(filenames=[], string='', scale=1, useDefaults=True, reload=False, simulation=0)
Initialize the XPLOR dihedral restraints (CDIH) potential term.
 
Parameters are:
   filenames   - either a single filename, or a sequence of filenames of
                 dihedral restraint assignment tables.
   string      - assignment table as a plain string.
   scale       - scale factor (defaults to 1).
   useDefaults - use the default sidechain restraints (default: True)
                 these force chi2 angles of PHE, TYR, ASP, and chi3 of GLU
                 to obey IUPAC naming conventions.
   reload      - if True, reload the previously-loaded dihedral restraints.
                 This is useful if the restraints are lost do to the
                 number of atoms changing (SCRATCHing).
initDynamics(ivm, bathTemp=-1, finalTime=0.2, numSteps=0, stepsize=0.001, potList=0, printInterval=50, initVelocities=0, eTol_factor=0.001, eTol_minimum=0)
Initialize an ivm.IVM object for PC6 (6th order predictor-corrector)
dynamics.
 
The default value of bathTemp is the ivm's current value.
 
In addition to the function arguments, the following IVM parameters are
initialized:
  maxDeltaE       -> 10000
  responseTime    -> 5
  adjustStepsize  -> True
  scaleVel        -> True
  resetCMInterval -> 10
 
Additionally, if initVelocities!=0, the initial velocities will be
randomized.
 
eTolerance is set by the following formula:
   eTolerance = eTol_factor * bathTemp + eTol_minimum
initHBDA(filename, scale=500, simulation=0)
Initialize the XPLOR hbda potential term
initHBDB(selection='not pseudo', restraintFile=None, prnfrq=0)
Initialize the XPLOR 
hbdb
potential term database hydrogen-bonding term - h-bonds are determined
dynamically. 
 
The selection argument specifies which atoms in included in the
calculation. Unfortunately, the current implementation allows only contiguous
ranges of residue numbers. So the selection specifies only the segment name
and the minimum and maximum residue number.
 
prnfrq can be changed from zero to periodically print out HBDB info during
dynamics and minimization.
 
By default, HBDB runs in ``free'' mode, in which hydrogen bonds are
determined on the fly. If you instead wish to use a fixed list of
h-bonds, specify a filename for the restraintFile argument. 
 
The HBDB term may have issues with hetereo multimers. Please use with
caution.
initMinimize(ivm, potList=0, printInterval=10, numSteps=500, maxCalls=20000, dEPred=0.001)
Initialize an ivm.IVM object for Powell minimization.
 
In addition to the function arguments, the following IVM parameters are
initialized:
  constrainLengths
  maxDeltaE
  eTolerance
  gTolerance
initNBond(cutnb=4.5, rcon=4.0, nbxmod=3, selStr='all', tolerance=0.5, repel=None, onlyCA=0, simulation=0, suppressException=False)
standard initialization of the XPLOR VDW term, configured to use the
repel potential. The XPLOR nonbonded potential term is described
here.
 
Note that cutnb should be greater than rmax + 2*tolerance, where rmax
is the largest vdw radius.
 
selStr specifies which atoms to use in nonbonded energy calculations.
 
If onlyCA is True, this string gets set to "name CA" - but this option is
deprecated.
 
The default value of repel is 0.8 for pre-3.1 version Xplor-NIH protein
and nucleic acid parameters, and 0.9 for more later versions. If a repel
value less than 0.9 specified to be used with post-3.0 parameter sets, an
exception will be raised, unless suppressExceptions is set to True.
initOrie(system='dna', pairs=[], selection='all', distCutoff=1000)
Set up the nucleobase-nucleobase positional knowledge-based potential.
 
The argument system is either 'dna' or 'rna'.  pairs is a sequence
of (x, y) tuples, where x and y are the residue numbers of two
bases involved in a Watson-Crick (WC) interaction.  Alternatively,
pairs may be the name of a file specifying the residue numbers in
WC pairs; the format of the file should be, e.g.,
 
1  24
2  23
 
where residue 1 is paired with 24, and 2 with 23. The use of pairs requires
that the specified residue numbers be unique. i.e. not replicated in
separate segids.
 
The selection argument can be used to specify a subset of the system to
apply this term to, but its selection string is not consulted if the
pairs argument is specified.
 
distCutoff is used to limit calculated interactions in the nonconsecutive 
contributions to the ORIE term.
initParams(files=[], string='', reset=0, use_dihe=False, weak_omega=False, simulation=0, system=None, verbose=True, silent=None)
files is a structure type or filename or a list of a combination of
these two.
 
valid structure types are: protein, nucleic, and water
 
string is an optional string of XPLOR-formatted parameter statements.
 
The default protein, nucleic and water parameter files are specified in
the protocol.parameters dictionary. Default values can be modified from
a script using e.g.
   import protocol
   protocol.parameters['protein'] = '/path/to/other/file.par' 
 
If the argument is a filename, the current directory is first searched
for the file, and then the TOPPAR directory is searched.
XPLOR parameters are documented here.
If reset is true, all previous parameter settings are discarded
If weak_omega is true, improper force constants associated with the
peptide bond are reduced by 1/2 to allow a small amount of flexibility. If
weak_omega is set to a float, these force constants are scaled by that 
value. 
 
If use_dihe is True, the DIHE XPLOR term is enabled.
 
If system is specified, specify that the parameters are forthe specified
system. Valid values are 'protein', 'nucleic', 'water', and 'metal'.
 
If verbose is False, XPLOR messages are suppressed. The deprecated 
parameter silent has the opposite meaning.
initPlanarity(filenames=[], string='', basePairs=[], simulation=0)
Initialize the XPLOR planarity restraint (PLAN) potential term.
 
Parameters are:
    filenames - either a single filename (a string), or a sequence of such
                filenames of planarity restraint tables.
    string    - table as a plain string.
    bairPairs - a sequence of pairs of residues for which to
                generate planarity restraints. Each base pair can be
                represented as a pair of residue numbers, or as a pair of
                tuples specifying (segid,resid).
 
The restraint table(s) should only contain the restraints-planar-statement.
Each statement is of the type:
 
    group
       selection=<selection>
       weight=<real>
    end
 
which adds a new atom group to the planar restraints database, where
<selection> is an XPLOR selection involving the atoms (at least four)
restrained to form a plane, and <real> is the restraint's weight.
 
Note that this function will clear any previously set up planarity
restraints.
initRamaDatabase(type='protein', simulation=0, selection='all')
Standard initialization of the RAMA torsion-angle database potential.
 
The argument type (a string) is either 'protein' or 'nucleic', specifying
the system type (protein or nucleic acid, respectively).
 
The selection argument (a string) can be used to specify a
subregion of the molecular system in the case that one does not
wish the RAMA database to apply to all residues. This option is
only implemented for nucleic acids.
initRandomSeed(seed=None)
set the initial random seed. If the seed argument is omitted, it is set
to seconds from the epoch mod 10^7
initStruct(files=0, string=None, erase=True, verbose=True, simulation=0)
Read XPLOR PSF files.
 
The files argument is a filename (a string) or sequence of filenames to read.
First the current directory is searched for the file, and then TOPPAR
is searched.  Alternatively, a psf entry (starting with PSF...) can be
directly passed as the files argument. In addition, a string value of the
PSF can be passed via the string argument.
 
Any pre-existing structure information is erased unless the erase argument
is set to False.
initTopology(files=[], reset=0, simulation=0, verbose=True, system=None, useDIHE=False)
file is a system or filename, or a list of a combination of
these two.
 
Valid systems are: 'protein', 'nucleic', 'water', and 'metal'.
 
The default protein, nucleic and water topology files are specified in
the protocol.topology dictionary. Default values can be modified from
a script using, e.g.
   import protocol
   protocol.topology['protein'] = '/path/to/other/file.top' 
or by specifying a file and system in this function's arguments.
 
First the current directory is searched for the file, and then TOPPAR
is searched.
XPLOR topology is documented here.
If reset is true, all previous topology settings are discarded.
 
A non-current simulation.Simulation can be specified with the 
simulation argument.
 
If system is specified, specify that the topology file is for the 
specified system. Valid values are 'protein', 'nucleic', 'water', and 
'metal'.
 
As a side-effect, atomSelLang.abbreviations are added for the value
[system], if it is specified in the files or system arguments.
 
If useDIHE is True, delete all instances of '!DIHE_' from the topology file.
This has the effect of enabling all dihedral terms.
initialRandomSeed()
return the initial random seed value, or that from
simulationWorld.random if initRandomSeed has not been called.
learnParams(selection='all', which='bond angle improper', bondForceConstant=1000, angleForceConstant=500, improperForceConstant=500, outFilename='new.par')
Read covalent parameters from the currently loaded structure. By default,
equilibrium values for bond, angles, and improper dihedral angles are
learned, but a subset of these can be specified using the which argument,
which is a space-separated string of values. This learn prodecure does not
determine force constants- they are specified explicitly using the
appropriate argument. This function requires that coorindates be read,
and that the associated covalent bond, angles, and impropers have been
defined.
loadPDB(file='', string='', entry='', model=-1, simulation=None, correctSymmetricSidechains=True, deleteUnknownAtoms=False, maxUnreadEntries=(20, 0.2), processSSBonds=True, processBiomt=False, allowedAltLoc=None, loadOccupanciesBfactors=False, verbose=-1, **kwargs)
Load a PDB entry from file, string, or from the wwPDB database.
autogenerate the PSF information and initialize coordinates.
 
deleteUnknownAtoms and maxUnreadEntries are passed to initCoords().
 
processSSBonds specifies whether the SSBOND pdb record is used in PDB
generation, and is passed to psfGen.pdbToPSF.
 
processBiomt specifies whether the REMARK 350 BIOMT record is used. By
default, Biomolecule 1 is read, specifying a different integer to the
processBiomt argument will read that entry.
 
loadOccupanciesBfactors is passed to initCoords.
 
The following optional arguments are passed through to other functions:
 
  swapProtons and assumeSegid are handled by the function
  matchInexactAtomEntry.
 
  includeHETATM, fixMethylImpropers and fixProtonImpropers, and 
  allowedAltLoc are handled by the function initCoords.
 
  useSeqres, deprotonateHIS, suppressExceptions and failIfInsertion are 
  passed through to psfGen.pdbToSeq.
 
  beginPatch, endPatch, convertToGly and useVariantResnames are passed 
  to psfGen.seqToPSF.
 
  deleteMissingResidues,  and deleteMissingAtoms are passed to 
  psfGen.pdbToPSF() to remove atoms as specified in PDB headers.
 
A side-effect is that protein and/or nucleic force field parameters are
loaded by initParams().
 
This function returns the InitCoordsResult from initCoords (described
below).
 
The verbose argument can be a bool or integer, with 0 or False resulting
in no output, and larger value producing more. The default value of -1
specifies that simulationWorld.world().logLevel() is consulted, and
that verbose is set to True if logLevel != "none".
massSetup(atomicMass=100, friction=10)
set uniform mass and friction values for all atoms.
Pseudo atoms may have nonuniform masses.
matchInexactAtomEntry(entry, pdbTool, selection, selectionSpecified, strictResNames, assumeSegid=None, swapProtons=False, verbose=0)
try to match a PDB record to an atom if the name is not exactly the
same. The coordinates will only be overwritten if they are unknown.
 
The current algorithm for matching is as follows:
  match H to HT1 
  else match if there's only a single undefined atom whose name starts
  with the same character.
  else
  match if atom names start with the same character and the remainder
  of the characters are simple tranpositions of each other.
  else
  match if there's only a single undefined atom whose name starts
  with and ends with the same character.
  else
  match if the first name of the name is a digit present in the correct
  atom name and the second character matches the first character of the
  correct name.
  else
  match if the two names have two matching characters.
  else
  match up   O --> OT1
           OXT --> OT2
 
  If the match is not unique, it is not made.
 
The residue number and segment id must match.
 
Entries with nonnumeric residue numbers are ignored.
 
Returns 0 1 or 2. If an atom is matched, 1 is returned, If an atom is
not matched, and a selection is not specified, then 0 is returned. If
a non-"all" selection is specified, or for entries whose residue names
correspond to those of pseudo atoms (e.g. ANI), a 2 is returned.
This function prints a warning on a non-matching entry if the selection
is not specified.
 
If strictResNames is True, then residue names in the PDB must match
those in the PSF (default), otherwise, the residue names need not match
and backbone atoms for mutants can be read. Note that sidechain atoms (or
base atoms in nucleic acids) may be strangely placed if this option is
False.
 
If the argument swapProtons==True, this function has the following
side effect: for gylcine, which should have HB1 and HB2, but
the input PDB names them HB2 and HB3, stereo-assignment is maintained by
swapping the positions. swapProtons defauls to False.
 
The argument assumeSegid affects matchInexactAtomEntry's behavior when
the segid field is blank. If the PSF has a blank segid, and the PDB has
assuemSegid, the atom will match. Likewise, if the PDB has a blank segid
(and chainID), the entry will match if the PSF has assumeSegid.
splitModel(filename, explicitModel=-1, explicitAltLoc=None)
Process a filename with a trailing model number or altLoc character
separated by a ``:'' character from the preceeding filename. The field
following the colon is interpreted as a model number if numeric, and
as a non-space altLoc value if alphabetical, in which case it should be
a single character. 
 
The value of explicitModel overrides the colon-separated model if it is 
not negative. The value of explicitAltLoc overrides a colon-separated 
altLoc if not None.
 
Returns a tuple of (file name,model, altLoc)
syncArrayBeforeDelete(anArray, selectionToDelete)
Given an array of values indexed by atom number, and a selection of
atoms to be deleted, remove the entries in the array so that they will
match up with the correct atoms after the atoms are actually deleted.
topologyFilename(system)
Function used to safely query the topology filename associated with the
specified system. It re-initializes the topology entry associated with
system to a TopologyEntry object if it is a bare string.
torsionTopology(ivm, fixedOmega=False, oTensors=[], breakRiboseSel="name C4' or name O4'", breakProlineSel='name CD or name CG', flexRiboseRing='[nucleic]')
Configure the ivm.IVM topology for standard torsion angle setup.
 
Given input ivm (an ivm.IVM instance), this function groups
rigid side chains and "breaks" disulfide bonds and proline and
ribose rings.
 
The bond broken in proline and ribose rings is respectively
specified by arguments breakRiboseSel and breakProlineSel, each an
XPLOR-style selection string selecting the atoms in the bond.  The
lengths of broken bonds are to be restrained by the force field as
they do not take on fixed values like those of "unbroken" bonds.
Bond angles and improper dihedrals associated with the atoms of
broken bonds are also flexible, and need to be restrained by the
force field.  In addition to torsion angles, all endocyclic bond angles
within ribose rings can be made active by specifying an atomSelLang
selection string in flexRiboseRing (default: empty string) using
selectTools.IVM_flexibilizeRiboses. Such selection should involve
at least one atom from the targeted residues, as in the following
examples:
 
'not pseudo and tag'  -> one atom per residue, excluding pseudoatoms
'tag and resid 1:33'  -> one atom per residue in the specified residue
                         range
'tag and segid "A"'   -> one atom per residue in the specified segid
 
(Selection of atoms within residues that do not have a ribose ring will
cause the program to fail.)
 
If fixedOmega is set to True, also fix protein omega backbone angles.
 
oTensors is a list of varTensor.VarTensor objects used for alignment
tensors.
 
Warning: This function sets all groupings and hinge types which are not
already specified, so it must be called last in the setup of an IVM's
topology.
updatePseudoAtoms(sim=None)
Update the set of pseudo atoms selected by the PSEUdo atom selection
keyword.
writeCIF(filename=None, selection='all', entryName='XplorNIH-Structure', modelCoords=None, occupancies=None, bFactors=None, remarks=[])
Write atoms from specified selection to an mmCIF file with the specified
filename.
 
The selection argument should be an atomSel.AtomSel object specifying
atoms to write out.
 
By default, no entries are created for occupancies or bFactors- these 
can be included in output .cif files by filling the appropriate arguments.
If they are instead set to the string "default", then the associated
field is filled with a ".".
 
The contents of the written mmCIF file are returned.
 
No backup file is created in this version.
writePDB(filename, writeChainID=False, selection='all', occupancies=None, bFactors=None, saveOccupanciesBfactors=False, remarks=[], writeConect=False)
write atoms from specified selection to a PDB file with the specified
filename
 
If writeChainID is specified, the segid is stored in the chainID column.
 
The selection argument should be an atomSel.AtomSel object specifying
atoms to write out.
 
The arguments occupancies and bFactors, if specified, should be numerical
sequences of length len(selection). These are placed in the occupancies and
bFactor PDB fields.
 
If saveOccupancyBfactors, occupancies and bfactors are read from the
associated xplorSimulation.XplorSimulation XPLOR process and placed in
the appropriate pdb columns.
 
remarks should be a sequence of strings without newlines, which are added
as REMARKS statement entries to the PDB file.
 
If writeConect is True, bonding information for all atoms in selection
will be written out in CONECT records.

 
Data
        defaultSelectionString = 'all'
downloadPDB_connectAttempts = 3
massFuncs = []
parVersion = {}
parameters = {'axis': 'axes.par', 'metal': 'ion.par', 'nucleic': 'nucleic.par', 'protein': 'protein.par', 'water': 'tip3p.parameter'}
parametersInitialized = {'axis': [], 'metal': [], 'nucleic': [], 'protein': [], 'water': []}
pdbLocation = 'https://files.wwpdb.org/pub/pdb/data/structures/divided/pdb/'
pseudoResNames = {'ANI', 'TAU', 'TEMP'}
systems = ('protein', 'nucleic', 'water', 'metal', 'axis')
topVersion = {}
topology = {'metal': TopologyEntry(filename='ion.top', linkages=None,...Statements=None, lastStatements=None, options={}), 'nucleic': TopologyEntry(filename='nucleic.top', linkages=N...Statements=None, lastStatements=None, options={}), 'protein': TopologyEntry(filename='protein.top', linkages=N...Statements=None, lastStatements=None, options={}), 'water': TopologyEntry(filename='tip3p.topology', linkage...Statements=None, lastStatements=None, options={})}
topologyFuncs = []
topologyInitialized = {'metal': [], 'nucleic': [], 'protein': [], 'water': []}