solnScatPotTools
index


 
tools to aid in setup/analysis of potential terms employing solution
scattering data.
 
this module provides functions to simplify the analysis of
solnScatPot.SolnScatPot potential terms.

 
Classes
       
builtins.object
Interp

 
class Interp(builtins.object)
    Interp(listOfTuples)
 
class for linear interpoation (/extrapolation)
 
  Methods defined here:
__call__(s, x)
return the result of interpolation (/extrapolation) to the point x.
__init__(s, listOfTuples)
specify a list of (x,y) pairs sorted such that the x values are
monotonically increasing.

Data descriptors defined here:
__dict__

 
dictionary for instance variables (if defined)
__weakref__

 
list of weak references to the object (if defined)

 
Functions
       
Pr(numBins, rMax, selection='known', weights=None)
return the unnormalized pairwise distance distribution function.
 
numBins specifies the granulaity of the distribution, and rMax specifies
the upper distance cutoff.
 
The calculation is performed using all atoms in selection
  [all atoms whose positions have been defined, by default].
 
The optional weights array specifies a weight for each atom.
SolventAtomParams(volume=-1, radius=-1)
analyze(potList)
perform analysis of SolnScatPot terms and return nicely formatted
summary.
A Chi^2 value is calculated based on the assumption that the weights are
1/sigma_i^2, where sigma_i is the error in the observed value of I_i.
 
The deviation value should be nonzero only for EnsembleSimulation
calculations with Ne>1. It is calculated as
 
Dev^2 = 1/Ne 1/Nk \sum_i \sum_j w_j * (I_ij-I_j)^2
 
where Nk is the number of datapoints and all values of I are normalized by
calcdScale, and multiplied by globCorrect.
analyze_chi2(term)
compute the Chi^2 value for the solnScatPot.SolnScatPot term
argument.
chi2(weights, I, expt, normalizeIndex=-3)
return chi^2 using the specified normalization method
fitParams(scat, r0Start=None, r0End=None, r0Num=None, V0Start=None, V0End=None, V0Num=None, rhobStart=None, rhobEnd=None, rhobNum=None, bgStart=None, bgEnd=None, bgNum=None, writer=None, verbose=False)
given a SolnScatPot instance, optimize the boundary layer
solvent density, the excluded volume effective radius, and excluded
volume parameters. These parameters are the same as those optimized in
the crysol program (J. Appl. Cryst. 28, 768-773 (1995)]
 
scat is a SolnScatPot instance
 
   r0Start=0.9,  r0End=1.1,    r0Num=21,
   V0Start=0.92, V0End=1.08,   V0Num=21,
   rhobStart=0,  rhobEnd=0.02, rhobNum=20,
   bgStart=-0.5, bgEnd=0.5,    bgNum=0,
 
specify the ranges for the three fitting parameters r0, V0 and rhob,
and a constant background present in the experimental spectrum. The r0
and V0 parameters specify a multiplicative factors relative to rm and
Vm, respectively. By default, bg=0 and no background fitting is performed.
bgStart and bgEnd specify values relative to the average value of the
last 10 experimental datapoints. Default values of all parameters are
specified above. 
 
If writer is specified, it should be a function to take informational
output.
fitSolventBuffer(scat, qValues=None, sampleData=None, bufferData=None, Iweights=None, minQ=None, maxQ=None, alpha0=1, sigmaAlpha=0.01, sigmaRho0=0.01, sigmaRhob=0.1, sigmaBG=0.01, rhob0=0.001, bg0=0, N0=None, optimizeRho0=False, optimizeRhob=True, writer=None, computeContributions=False, verbose=False)
given a SolnScatPot instance, and input sample and buffer filenames,
optimize the buffer subtraction, background contribution and boundary
layer contribution
 
scat         -a SolnScatPot instance.
sampleData   -a sequence containing the sample scattering curve, without
              buffer subtraction. The length must be the same as
              scat.expt().
bufferData   -a sequence containing the buffer scattering curve. The
              length must be the same as scat.expt().
 
This is an adaptation of the AXES algorithm:  A. Grishaev, L.A. Guo,
T. Irving, and A. Bax, ``Improved fitting of solution X-ray
scattering data to macromolecular structures and structural
ensembles by explicit water modeling'',  J. Am. Chem. Soc. 132,
15484-15486 (2010).
 
 
The experimental curve is computed as:
 
   I_expt(q) = I_samp(q) - alpha I_buff(q) + bg
 
The calculated curve is given by:
 
   I(q) = N * < |A(q) + rhob Ab(q)|^2 >
 
where N is overall normalization, A(q) is the scattering amplitude from
solute molecule and excluded solvent,  Ab(q) is scattering from
the solvent boundary layer, and rhob is the coefficient of the boundary
layer scattering.
 
The fit range can be reduced from the entire range by specifying
minQ or maxQ.
 
Initial values of the following parameters may be specified: 
 
    alpha0     - initial value of coefficient for buffer subtraction
    sigmaAlpha - uncertainty in alpha
    sigmaRho0  - uncertainty in rho0 - only used when fitting rho0.
    bg0        - initial value for isotropic background
    N0         - initial value of normalization for the calculated curve.
                 By default, this is taken by an overall chi^2 fit to
                 the initial experimental curve.
 
The density of excluded solvent can be optimized by setting
optimizeRho0=True.
 
The computation of the optimal bound-solvent density can be disabled by
setting optimizeRhob=False.
 
The return value is a dictionary with optimal values in elements
indexed by 'alpha', 'bg', 'N', 'rho0', and 'rhob', in addition to
'chi2' and 'cost', the total function value where the parameters
take their optimal values. If the computeContributions argument is
set to True, the dictionary will additionally include 'Ivac',
'Isol', 'Ib'and 'Inob', corresponding to intensities from the
vacuum, excluded solvent boundary, and boundaryless contributions.
 
If writer is specified, it should be a function to take informational
output.
interpolateCurve(experiment, weights, numPoints)
given a curve with elements (q,I) and corresponding weights, return the
corresponding linearly interpolated curves sampled uniformly over
numPoints.
normalize(weights, I, expt, normalizeIndex)
readData(experiment, preweighted=False, minQ=None, maxQ=None)
read SAXS or SANS data from file, string, or sequence. If a
sequence, each element should be a sequence of 2 or more values;
the 1st two are q,I, while the 3'rd is an optional weight
value. If file, it is read and processed as a string. If a string,
leading # or ! characters start comment lines which are ignored.
Noncomment lines must contain at least 2 space-separated values
corresponding to q and I. An optional third value specifies experimental
error or weight.
 
The preweighted and maxQ args specify how string (and file) data is
interpreted. preweighted=True specifies that the 3rd value is a weight,
while a False value specifies that it is an experimental error (and
a weight of 1/e^2 is computed). If maxQ is specified, don't read data with
larger q values. If minQ is specified, don't read data with
smaller q values.
 
Returns a list of (q,I) tuples and a list of weights for each of these.
reduce(...)
reduce(function, iterable[, initial]) -> value
 
Apply a function of two arguments cumulatively to the items of a sequence
or iterable, from left to right, so as to reduce the iterable to a single
value.  For example, reduce(lambda x, y: x+y, [1, 2, 3, 4, 5]) calculates
((((1+2)+3)+4)+5).  If initial is present, it is placed before the items
of the iterable in the calculation, and serves as a default when the
iterable is empty.
splineCurve(experiment, weights, numPoints)
given a curve with elements (q,I) and corresponding weights, return the
corresponding splined curves sampled uniformly over numPoints.
useGlobs(term, globTable=[], globRules=[], verbose=0, weightFunction=<function <lambda> at 0x7f2aff557d80>)
set up solnScatPot.SolnScatPot term to use the atom globbing
approximation.
 
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 atom names to be globed.
 
weightFunction takes an atom as an argument and returns the relative
weight to give it within the glob.
 
Atoms in term.selection() which are not specified by globTable or by
globRules are placed into single-atom globs.

 
Data
        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')], ...}
params = {'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), ...}
solventParamSets = {'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), ...}}
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), ...}}