solnScatPot |
index |
potential term for solution-phase scattering experiments, including
SAXS and SANS
One normally creates SolnScatPot objects using
solnXRayPotTools.create_XRayScatPot or sansPotTools.create_SANSPot.
This term is ensembleSimulation.EnsembleSimulation - aware.
See C.D. Schwieters and G.M. Clore, ``Using Small Angle Solution
Scattering Data in Xplor-NIH Structure Calculations ,''
Progr. NMR Spectroscopy, accepted (2014).
constructor:
SolnScatPot(instanceName ,
sel ,
qValues ,
groupMap ,
formFactors ,
iFormFactors ,
solventParams,
rho0 ,
radiusScale ,
volumeScale ,
globDefs ,
experiment ,
qValuesExpt ,
simulation )
where:
instanceName user-specified identifier.
sel atomSel.AtomSel which specifies atoms to use in
the scattering calculation (frequently all heavy atoms).
qValues a sequence of values of scattering vector at which
to perform refinement. q=4 PI sin(theta) / lambda,
where 2 theta is the scattering angle, and lambda is
the wavelength of the incident beam. These should usually
be specified at uniform values of q.
groupMap a sequence of group names to use in the scattering calc.
formFactors a sequence of maps of group name to effective form factor,
one for each value of scattering vector- in the same order
iFormFactors a sequence of maps of atom index to effective form factor,
one for each value of scattering vector- in the same order.
Generally, individual atoms need not be specified-
but are rather specified by group name. This is for
exceptions, and will override any group entry.
solventParams a dictionary of 2-membered tuples- one for each atom group
name. Each tuple corresponds to the radius and volume for
that atom group.
rho0 excluded solvent electron density.
radiusScale scale factor used in the calculation of the displaced
solvent scattering contribution. See below.
volumeScale scale factor used in the calculation of the displaced
solvent scattering contribution. See below.
globDefs a sequence of glob definitions. Each definition is a
sequence of tuples (w,index), where w is a weight and index
is an atom index.
qValuesExpt The q-values corresponding to the experimental measurements.
The fit of calculated curve to experiment is done for
each value of q here, with the calculated curve interpolated
using a cubic spline, hence it is important that none of
the qValuesExpt fall outside the range given in qValues.
experiment solvent-corrected experimental scattering intensity values-
one for each value of q in qValuesExpt.
simulation optional simulation.Simulation object.
members:
The following parameters can be set [defaults in square brackets]
scale - scale factor (force constant) [1]
weights - sequence of factors to weight calculated scattering intensity
to its corresponding experimental value. Used for rms, energy
calculation. See below. Default is a value of 1 for all.
calcType - specifies how I(q) is to be calculated. Valid values are
"n2" - use the exact, but slow Debye formula.
"n" - numerically integrate over surface of constant
q, using theta, phi values determined by angleType
"uniform" - a faster version of the n algorithm which
can be used if q values are equally spaced.
"multipole" - Use the spherical harmonic expansion of
Svergun, Barberato and Koch,
J. Appl. Cryt. 28, 768 (1995)
Note: at this time it is not possible
to refine using the calcType, as the
derivatives are not yet coded.
numAngles - for calcTypes n and uniform, the number of (theta,phi)
pairs to use in the evaluation of I(q).
qValues - the values of scattering vector amplitude q.
expt() - the observed values of I(q).
formFactors - specified the formFactor map.
iFormFactors - specify the atom-indexed formFactor map.
angleType - specifies how (theta,phi) angle pairs are chosen over
the surface of a sphere. Choices are:
"random" - choose randomly- not recommended.
"fibonacci" - method based on Fibonacci numbers.
D. Svergun, Acta Cryst. A50, 391 (1994).
"spiral" - spiral algorithm from
E.B. Saff and A.B.J. Kuijlaars,
Distributing Many Points on a Sphere,
The Mathematical Intelligencer, 19(1),
Winter (1997).
This is the preferred algorithm.
lMax - for calcType multipole, the maximum value of l to use in
the spherical harmonic expansion.
cmpType - should be "plain" or "log." This controls whether the
logarithm is taken before the energy is evaluated. See below.
[Default: plain]
globs - sequence of glob definitions - see description above.
normalizeIndex - if greater than or equal to zero, specifies which
value of q at which to consider Icalc and Iobs to be
equal. If the value is -2, the normalization is taken
as the average of the Icalc or Iobs values weighted by
weights. If the value of normalizeIndex lies outside
the above values, no normalization is performed. If the
value -3 is specified, the normalization will instead
minimize chi^2.
[default: -3]
radiusScale - scale factor used in the calculation of the displaced
solvent scattering contribution. See below.
volumeScale - scale factor used in the calculation of the displaced
solvent scattering contribution. See below.
radiusType - "radius" or "volume" - specifies the method used for
computing atomic radii. See below.
rho0 - excluded solvent electron density.
calcBoundary - boolean - whether to calculate boundary layer
scattering contribution [default: False]
rhob - boundary layer electron density [default: 0.03]
boundaryThickness - the thickness of the boundary layer [default: 2.88].
solventRadius - radius used in the surface computation
[default: 1.44].
bg - isotropic background value to be subtracted from
observed scattering values.
formScale - overall scale factor for vacuum atomic form factors.
useGlobs - whether or not to use globs in the scattering calculation.
verbose - if true, produce more detailed output [False]
useSimEnsWeight - whether to use the ensemble wieghts set with setEnsWeights
or to use those of the underlying EnsembleSimulation.
ensWeights - a sequence of ensemble weights to use when calculating
the ensemble- averaged value of each scattering curve. By
default, the weights of the underlying EnsembleSimulation
are used. If these are overridden by calling the
setEnsWeights method, and useSimEnsWeight is set to False.
the above quantities may be retrieved using the member function form
quantity(), while they are set using the form setQuantity(value).
QValues - a plain member containing values of (x,y,z) on a r=1
sphere used for numerical integration of I(q).
methods:
calcEnergy() - calc energy, returns the energy value.
calcEnergyAndDerivs(derivs) - calc energy, derivs, returns the energy value.
calcGlobCorrect(val) - calculate the globic correction factor. The argument
val can be calcType (e.g. 'n2'), in which case
the globic correction factor will be calculated by
using the specified calcType. val can alternatively
be a raw array of I(q) values which are used to
generate g(q). val defaults to the current calcType.
globCorrect() - return the current globic correction factor.
getF() - return a matrix of scattering amplitudes at each angle
and q value.
rms() - return the rms of calcedShift - effShift
info() - string with current info about the state of this instance
simulation() - return the associated simulation.Simulation.
calcd() - return the calculated values of I(q) for each value of
qValues.
splined() - return the calculated values of I(q) for each value of
qValuesExpt.
exptScale() Iobs(qn), for n=normalizeIndex, if normalizeIndex>-1,
1 otherwise.
calcScale() Icalc(qn), for n=normalizeIndex, if normalizeIndex>-1,
1 otherwise.
I_contrib(memberIndex) - the scatter intensity from the specified ensemble
member structure, or the local member if the
index is omitted.
boundaryI() - scattering intensity contribution from the boundary
layer.
selection() - return the atom selection specified in the
constructor.
radii() - return array of atomic radii
aveRadius() - average atomic radius.
aveVolume() - average atomic volume.
boundaryVol() - volume of the boundary layer.
calcBoundaryPoints() - calculate new points used to compute the
boundary layer scattering contribution.
ensWeight(index) - return the ensemble weight associated with the
specified member.
addRigidRegion(atomSel) - indicate that the specified atoms move as
a rigid body so that a very fast update
algorithm can be used. Use of this method
does not actually constrain the relative
atom positions to be constrained. For
that, the atoms must be grouped in the
appropriate ivm.IVM objects.
solnScat() - return a SolnScat object. This object is
used to calculate the raw scattering curve
from atomic coordinates. It has a few
useful methods to access properties not
directly accessible from the SolnScatPot
class:
numRigidRegions() - return the number of rigid regions.
rigidRegion( n ) - return the definition of rigid region
n, indexed from 0.
Accessor members:
innerBoundaryFilename - if set, write out the inner molecular surface
used for boundary contribution calculation.
outerBoundaryFilename - if set, write out the outer molecular surface
used for boundary contribution calculation.
Energy is calculated as
scale * rms^2
where
rms^2 = sum_i w_i (Icalc_i - Iobs_i)^2 / (Nq),
in which Nq is the total number of values of q at which Iobs is
specified, and w_i is the corresponding weight. Typically, we choose
w_i = (sigma_i)^(-2), where sigma_i is the error in I(q) at q=q_i.
Given a set of atomic coordinates the scattering intensity is
calculated from the following expression:
I(q) = <| A(Q) |^2>_q
where <>_q denotes average over solid angle at constant amplitude of
the scattering vector.
The scattering amplitude A(Q) consists of a sum of contributions from
each atom (usually only heavy atoms are included)
A(Q) = sum_j f_j(q) exp(i dot(Q,r_j)) +
sum_j2 rhob(r_j2) exp(i dot(Q,r_j2))
where f_j(q) is the total effective scattering amplitude of atom j
(including excluded solvent effects), and r_j is atom j's position. The
second contribution is due to boundary-layer scattering, and the sum
over j2 is, in principle an integral over the region in which rhob,
the boundary layer density is nonzero.
For rhob=0, the integral over solid angle can be carried out
analytically, resulting in the Debye formula calcType=n2:
I(q) = sum_j1 sum_j2 f_j1(q) f_j2(q) sinc(q |r_j1-r_j2|)
where sinc(x) = sin(x) / x.
This expression is quite expensive for biological systems due to the
double sum.
One alternative is to express A(Q) as an expansion in terms of
spherical harmonic functions (calcType=multipole), and approximate
I(q) by truncating the expansion. This is works well (is efficient)
for small scattering angles - like those used in most SAXS
experiments. For larger angles, the expansion fails spectacularly.
One can also evaluate the integral over solid angle by simple grid
integration (calcTypes=n, uniform). This works well for all moderate
values of q.
Excluded Solvent Treatment:
The total atomic effective scattering amplitude is
f_j(q) = fv_j(q) - g_j(q),
where fv_j(q) is the vacuum atomic form factor, and g_j is an excluded
solvent contribution
g_j(q) = rho0 * V_j * exp(-s**2 * PI * V_j**(2./3)) *
volumeScale * exp(-q**2 * PI * (4PI/3)^{1/3} * (r0^2 - rm^2)),
which depends on atomic volume V_j and scaled scattering vector
amplitude s=q/(2PI). rm is aveRadius() and r0 is a scaled radius
r0= radiusScale() * rm
If radiusType=="volume" (the default), rm actually is computed as the
radius corresoponding to aveVolume() [i.e. (3 Vm / (4 PI))^(1/3) ]
Boundary Layer Treatment:
The boundary layer contribution is computed based on a molecular
surface generated using surfD. The molecular surface is generated
based on atomic radii, boundaryThickness and solventRadius(). An outer
surface is generated by rolling a ball of radius solventRadius over
the surface generated of spheres of radius
radii+boundaryThickness. Then, an inner surface is generated by, for
each surface triangle, dropping a line segment of length
boundaryThickness in the direction opposite the surface normal. Each
surface triangle thus generates a voxel. The contribution of each
voxel to the scattering amplitude is represented by a sphere sharing
the center and volume of the voxel.
Ensemble Calculations
When calculated in the context of an ensembleSimulation.EnsembleSimulation
the scattering intensity becomes:
I(q) = sum_n W_n I_n(q)
where the sum is over all members of an ensemble, w_n is the weight
for ensemble member n, and I_n(q) is the scattering intensity from
ensemble member n, calculated as above.
One should note that the approximate methods of calculating I(q)
introduce a rotational dependence: the grid on which the scattering
amplitude is calculated is fixed in space. This rotational dependence
implies that the potential energy term is not constant upon overall
rotation. The amount of the rotational dependence decreases with
increased numAngles. A future enhancement may have the grid rotate
along with the rest of the system.
Helper functions:
The following are used in solnScatPotTools:
pairDistribution(atomSel,
weights,
numBins,
rMax ) - compute the pair distance distribution due to
the positions of the atoms in the selection
atomSel, with each position weighted by
the corresponding element in weight (should
have length equal to number of atom in
atomSel). The returned histogram has numBins
bins with distance ranging between 0 and rMax.
volContrib(r0,
rm,
qValues,
Fs) - generate matrix
exp( c * qValues^2 ) * Fs
where c = -(4*PI/3)^(2/3)*(r0^2-rm^2)/(4*PI)
This is OMP-parallelized, and likely not safe
to call in threaded region.
calcIfromF(F,
pot) - return a splined value of an ensemble-
averaged I computed from F, where I for
ensemble member i is computed as the square
magnitude of the i-th row of F. The pot
argument is a SolnScatPot instance.
calc_dIdrhob(Fb, - return the derivative of I with respect to
F, rhob.
rhob,
pot)
# This file was automatically generated by SWIG (http://www.swig.org).
# Version 4.0.2
#
# Do not make changes to this file unless you know what you are doing--modify
# the SWIG interface file instead.
Classes | ||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
|
Functions | ||
|
Data | ||
ANGLE_FIBONACCI = 1 ANGLE_RANDOM = 0 ANGLE_SPIRAL = 2 __CDSVector_hh__ = 1 |