solnXRayPotTools
index


 
tools to aid in setup/analysis of potential terms employing solution
X-Ray scattering data.
 
this module provides functions to simplify the creation and
analysis of solnScatPot.SolnScatPot potential terms. 
 
See C.D. Schwieters and G.M. Clore, ``Using Small Angle Solution
          Scattering Data in Xplor-NIH Structure Calculations ,''
          Progr. NMR Spectroscopy, accepted (2014).

 
Functions
       
calcFormFactors(q, atomGroups)
calculate atomic form factors for combined atoms in atomGroups.
calcFormFactors_indexed(q)
calculate atomic form factors for special combined atoms identified by
an atom selection and references by atom index.
calc_g_jbw(rho0, volume_scale, V, q, r0, rm)
convertFormFactors_jbw(sel, fValues, qValues, groupMap, rho0, radiusType, radius_scale, solvent_volume, volum_scale)
create_SolnXRayPot(instanceName, aSelection='not PSEUDO and not name H*', experiment=0, radiusScale=1, volumeScale=None, radiusType='volume', rho0=0.334, normalizeIndex=-3, preweighted=True, numPoints=None, sampleData=None, bufferData=None, useInternalSpline=False, minQ=None, maxQ=None, formFactors='original', useNM=False, verbose=None, simulation=0)
create a solnScatPot.SolnScatPot term appropriate for refining
against solution X-Ray scattering data.
 
experiment can be specified as the name of a file whose first two
columns are the experimental q, I(q). An optional third column specifies
a weight to use for I(q). Lines starting with # or ! are
skipped. 
 
If preweighted=False, the third entry is standard deviation instead of
a weight, and the weight is computed as 1/sigma^2 if sigma!=0, else 0.
If normalizeIndex>=0, the weight is also normalized [if this behavior is
not desired, set normalizeIndex after this function has been called.]
 
experiment can also be specified as a sequence of (q,I(q)) tuples.
 
experimental value of q less than minQ will be ignored.
experimental value of q greater than maxQ will be ignored.
 
normalizeIndex specifies which grid point to use to normalize data.
normalizeIndex=-2 specifies average normalization .
normalizeIndex=-3 specifies normalization which minimizes the Chi^2 value.
Other negative values mean no normalization.
 
numPoints specifies the number of datapoints to sample. The experimental
I(q) and errors/weights are linearly interpolated and sampled uniformly
over the entire range. No smoothing is performed.
 
If the q values are evenly spaced or numPoints is specified,
calcType is set to 'uniform'.
 
solnScatPot.SolnScatPot.numAngles is initialized to 500
 
solnScatPot.SolnScatPot.cmpType is initialized to 'plain'.
 
radiusScale is used to scale the effective atomic radii- perhaps due to
the presence of a solvent layer.
 
If specified, volumeScale is a separate scale factor for the effective
volume. If it is not specified, this scale factor will be consistent with
radiusScale
 
rho0 is the solvent electron density. The default value is that of water.
 
this function initializes effective atomic form factors including
an excluded volume term.
 
For atom j
  feff_j(q) = f_j(q) - g_j(q)
 
where f_j(q) is the atomic scattering factor of group j, evaluated using a
5-Gaussian sum.
 
g_j(q) = rho0 * V_j * exp(-s**2 * PI * V_j**(2./3)) *
         volumeScale * exp(-s**2 * PI * (4PI/3)^{1/3} * (r0^2 - rm^2))
 
where s=q/(2*PI), rho0 is the solvent electron density, V_j is the
solvent-displaced volume. The second line of the expression is an overall
expansion factor, in which 
 
rm is the radius corresponding to the average atomic group volume
(radiusType=volume) or the average atomic group radius (radiusType=radius),
and r0 is an effective radius.
 
  r0 = radiusScale * rm
 
and radiusScale is specified as a number approximately equal to 1.
 
Specify useNM=True to indicate that q is in units of inverse nm rather than
the default inverse angstroms.
create_solnXRayPot = create_SolnXRayPot(instanceName, aSelection='not PSEUDO and not name H*', experiment=0, radiusScale=1, volumeScale=None, radiusType='volume', rho0=0.334, normalizeIndex=-3, preweighted=True, numPoints=None, sampleData=None, bufferData=None, useInternalSpline=False, minQ=None, maxQ=None, formFactors='original', useNM=False, verbose=None, simulation=0)
create a solnScatPot.SolnScatPot term appropriate for refining
against solution X-Ray scattering data.
 
experiment can be specified as the name of a file whose first two
columns are the experimental q, I(q). An optional third column specifies
a weight to use for I(q). Lines starting with # or ! are
skipped. 
 
If preweighted=False, the third entry is standard deviation instead of
a weight, and the weight is computed as 1/sigma^2 if sigma!=0, else 0.
If normalizeIndex>=0, the weight is also normalized [if this behavior is
not desired, set normalizeIndex after this function has been called.]
 
experiment can also be specified as a sequence of (q,I(q)) tuples.
 
experimental value of q less than minQ will be ignored.
experimental value of q greater than maxQ will be ignored.
 
normalizeIndex specifies which grid point to use to normalize data.
normalizeIndex=-2 specifies average normalization .
normalizeIndex=-3 specifies normalization which minimizes the Chi^2 value.
Other negative values mean no normalization.
 
numPoints specifies the number of datapoints to sample. The experimental
I(q) and errors/weights are linearly interpolated and sampled uniformly
over the entire range. No smoothing is performed.
 
If the q values are evenly spaced or numPoints is specified,
calcType is set to 'uniform'.
 
solnScatPot.SolnScatPot.numAngles is initialized to 500
 
solnScatPot.SolnScatPot.cmpType is initialized to 'plain'.
 
radiusScale is used to scale the effective atomic radii- perhaps due to
the presence of a solvent layer.
 
If specified, volumeScale is a separate scale factor for the effective
volume. If it is not specified, this scale factor will be consistent with
radiusScale
 
rho0 is the solvent electron density. The default value is that of water.
 
this function initializes effective atomic form factors including
an excluded volume term.
 
For atom j
  feff_j(q) = f_j(q) - g_j(q)
 
where f_j(q) is the atomic scattering factor of group j, evaluated using a
5-Gaussian sum.
 
g_j(q) = rho0 * V_j * exp(-s**2 * PI * V_j**(2./3)) *
         volumeScale * exp(-s**2 * PI * (4PI/3)^{1/3} * (r0^2 - rm^2))
 
where s=q/(2*PI), rho0 is the solvent electron density, V_j is the
solvent-displaced volume. The second line of the expression is an overall
expansion factor, in which 
 
rm is the radius corresponding to the average atomic group volume
(radiusType=volume) or the average atomic group radius (radiusType=radius),
and r0 is an effective radius.
 
  r0 = radiusScale * rm
 
and radiusScale is specified as a number approximately equal to 1.
 
Specify useNM=True to indicate that q is in units of inverse nm rather than
the default inverse angstroms.
exp(x, /)
Return e raised to the power of x.
f_glob(q, atoms, dist)
determine globbed value of f (using calc_f) given a list of atom names
and a dictionary containing atom-atom distances. Can only be used with very
simple glob defs.
genFormFactors(qValues, aSelection)
getAveRadius_jbw(sel, groupMap, radius_type, radius_scale, solvent_volume)
getCombinedAtoms2(atom, bondsByAtomIndex, atomGroups=0)
given an atom object, determine the combination of heavy and light
atoms in this group.
 
atomGroups is used to hold per-group counts.
numElectrons(atom)
return the number of electrons of a given atom.
sinc(x)
updateQValues(pot, newQ)
update the solnScatPot.SolnScatPot term with new values of scattering
vector amplitude newQ.
use4GaussianFormFactors()
use5GaussianFormFactors()
useGlobs(term, globTable=[], globRules={'GLY': [('C', 'O', 'N', 'HN'), ('CA', 'HA1', 'HA2')], 'ALA': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'HB3')], 'VAL': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB', 'CG1', 'HG11', 'HG12', 'HG13', 'CG2', 'HG21', 'HG22', 'HG23')], 'LEU': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2'), ('CG', 'HG', 'CD1', 'HD11', 'HD12', 'HD13', 'CD2', 'HD21', 'HD22', 'HD23')], 'ILE': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB', 'CG2', 'HG21', 'HG22', 'HG23'), ('CG1', 'HG11', 'HG12', 'CD1', 'HD11', 'HD12', 'HD13')], 'PHE': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2'), ('CG', 'CD1', 'HD1', 'CD2', 'HD2', 'CE1', 'HE1', 'CE2', 'HE2', 'CZ', 'HZ')], 'TYR': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2'), ('CG', 'CD1', 'HD1', 'CD2', 'HD2', 'CE1', 'HE1', 'CE2', 'HE2', 'CZ', 'OH', 'HH')], 'TRP': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2'), ('CG', 'CD1', 'HD1', 'CD2', 'NE1', 'HE1', 'CE2', 'CE3', 'HE3', 'CZ2', 'HZ2', 'CZ3', 'HZ3', 'CH2', 'HH2')], 'ASP': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2'), ('CG', 'OD1', 'OD2')], 'GLU': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'CG', 'HG1', 'HG2'), ('CD', 'OE1', 'OE2')], 'SER': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'OG', 'HG')], 'THR': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB', 'OG1', 'HG1', 'CG2', 'HG21', 'HG22', 'HG23')], 'ASN': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2'), ('CG', 'OD1', 'ND2', 'HD21', 'HD22')], 'GLN': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'CG', 'HG1', 'HG2'), ('CD', 'OE1', 'NE2', 'HE21', 'HE22')], 'LYS': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'CG', 'HG1', 'HG2'), ('CD', 'HD1', 'HD2', 'CE', 'HE1', 'HE2', 'NZ', 'HZ1', 'HZ2', 'HZ3')], 'ARG': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'CG', 'HG1', 'HG2'), ('CD', 'HD1', 'HD2', 'NE', 'HE', 'CZ', 'NH1', 'NH2', 'HH11', 'HH12', 'HH21', 'HH22')], 'HIS': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2'), ('CG', 'ND1', 'HD1', 'CD2', 'HD2', 'CE1', 'HE1', 'NE2', 'HE2')], 'MET': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'CG', 'HG1', 'HG2'), ('SD', 'CE', 'HE1', 'HE2', 'HE3')], 'CYS': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'SG', 'HG')], 'PRO': [('C', 'O', 'N'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'CG', 'HG1', 'HG2', 'CD', 'HD1', 'HD2')], 'CYT': [('P', 'O1P', 'O2P'), ("O5'", "C5'"), ("O3'", "C3'"), ("O4'", "C4'"), ("C1'", "C2'"), ('N1', 'C2', 'O2'), ('N3', 'C4', 'N4'), ('C5', 'C6')], 'GUA': [('P', 'O1P', 'O2P'), ("O5'", "C5'"), ("O3'", "C3'"), ("O4'", "C4'"), ("C1'", "C2'"), ('N1', 'C2', 'N2'), ('N3', 'C4'), ('C5', 'C6', 'O6'), ('N7', 'C8', 'N9')], 'ADE': [('P', 'O1P', 'O2P'), ("O5'", "C5'"), ("O3'", "C3'"), ("O4'", "C4'"), ("C1'", "C2'"), ('N1', 'C2'), ('N3', 'C4'), ('C5', 'C6', 'N6'), ('N7', 'C8', 'N9')], 'THY': [('P', 'O1P', 'O2P'), ("O5'", "C5'"), ("O3'", "C3'"), ("O4'", "C4'"), ("C1'", "C2'"), ('N1', 'C2', 'O2'), ('N3', 'C4', 'O4'), ('C5', 'CM', 'C6')]}, verbose=0)
set up solnScatPot.SolnScatPot term to use the atom globbing
approximation for solution X-ray scattering.
 
globTable contains a list of list of atoms with in user-defined globs.
Atoms not specified in globTable are glob'ed by the pre-residue definitions
in the globRules dictionary.
 
globRules is a dictionary whose keys are upper case residue names
each entry containing a list of list of atoms to be globed.
 
Atoms in term.selection() which are not specified by globTable or by
globRules are placed into single-atom globs.

 
Data
        bondsByAtomIndex = {}
corr_facts = [(0.0, 1.0), (0.2, 1.0), (0.22, 1.06), (0.24, 1.1042), (0.26, 1.1333), (0.38, 0.9815), (0.68, 1.0088), (0.78, 1.0092), (0.92, 1.1406), (0.94, 1.1524), (0.96, 1.1571), (0.98, 1.1548), (1.3, 1.1637), (1.32, 1.1837), (1.34, 1.209), (1.36, 1.2233), (3.0, 1.2375)]
globRules = {'ADE': [('P', 'O1P', 'O2P'), ("O5'", "C5'"), ("O3'", "C3'"), ("O4'", "C4'"), ("C1'", "C2'"), ('N1', 'C2'), ('N3', 'C4'), ('C5', 'C6', 'N6'), ('N7', 'C8', 'N9')], 'ALA': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'HB3')], 'ARG': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'CG', 'HG1', 'HG2'), ('CD', 'HD1', 'HD2', 'NE', 'HE', 'CZ', 'NH1', 'NH2', 'HH11', 'HH12', 'HH21', 'HH22')], 'ASN': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2'), ('CG', 'OD1', 'ND2', 'HD21', 'HD22')], 'ASP': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2'), ('CG', 'OD1', 'OD2')], 'CYS': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'SG', 'HG')], 'CYT': [('P', 'O1P', 'O2P'), ("O5'", "C5'"), ("O3'", "C3'"), ("O4'", "C4'"), ("C1'", "C2'"), ('N1', 'C2', 'O2'), ('N3', 'C4', 'N4'), ('C5', 'C6')], 'GLN': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'CG', 'HG1', 'HG2'), ('CD', 'OE1', 'NE2', 'HE21', 'HE22')], 'GLU': [('C', 'O', 'N', 'HN'), ('CA', 'HA', 'CB', 'HB1', 'HB2', 'CG', 'HG1', 'HG2'), ('CD', 'OE1', 'OE2')], 'GLY': [('C', 'O', 'N', 'HN'), ('CA', 'HA1', 'HA2')], ...}
metalAtoms = ['ZN', 'MN', 'MG']
pi = 3.141592653589793
raw_corr_facts = [(0.0, 0.9882), (0.02, 0.9908), (0.04, 1.0087), (0.06, 1.0226), (0.08, 1.0242), (0.1, 1.0161), (0.12, 1.0087), (0.14, 1.0092), (0.16, 1.0074), (0.18, 1.0115), (0.2, 1.0256), (0.22, 1.06), (0.24, 1.1042), (0.26, 1.1333), (0.28, 1.1394), (0.3, 1.1171), (0.32, 1.0862), (0.34, 1.0438), (0.36, 1.0128), (0.38, 0.9815), ...]
solventVolume = {'C': (16.52189547022142, 1.58), 'CH': (21.688370252755973, 1.73), 'CH2': (26.521848780380633, 1.85), 'CH3': (32.02486353433834, 1.97), 'Ca': (32.02486353433834, 1.97), 'Cu': (8.78452975554819, 1.28), 'Fe': (7.986447935410647, 1.24), 'H': (5.131448118842135, 1.07), 'H2O': (19.3, 1.5133666759556506), 'MN': (9.2027720799157, 1.3), ...}
solventVolumeSets = {'svergun': {'C': (16.52189547022142, 1.58), 'CH': (21.688370252755973, 1.73), 'CH2': (26.521848780380633, 1.85), 'CH3': (32.02486353433834, 1.97), 'Ca': (32.02486353433834, 1.97), 'Cu': (8.78452975554819, 1.28), 'Fe': (7.986447935410647, 1.24), 'H': (5.131448118842135, 1.07), 'H2O': (19.3, 1.5133666759556506), 'MN': (9.2027720799157, 1.3), ...}, 'tiede': {'Br': (26.52, 1.6824786822618245), 'C': (9.0, 1.1735616258720787), 'CH': (20.0, 1.5314461446813734), 'CH2': (21.0, 1.556556282112301), 'CH3': (33.0, 1.8096574578554767), 'Ca': (31.0, 1.7723342420376789), 'Cl': (22.45, 1.5915874724077386), 'Cu(2)': (9.2, 1.1821910653824932), 'Fe(2)': (8.3, 1.1423111811746613), 'Fe(3)': (8.3, 1.1423111811746613), ...}, 'xiaobing': {'Br': (26.52, 1.6824786822618245), 'C': (16.44, 1.4345829602682996), 'CH': (21.59, 1.5709991296307175), 'CH2': (26.74, 1.6871182809580803), 'CH3': (31.89, 1.7891355151247366), 'Ca': (31.89, 1.7891355151247366), 'Cl': (22.45, 1.5915874724077386), 'Cu(2)': (9.2, 1.1821910653824932), 'Fe(2)': (8.3, 1.1423111811746613), 'Fe(3)': (8.3, 1.1423111811746613), ...}}