5. Small angle scattering (sas)¶
This module allows smearing/desmearing of SAXS/SANS data.
Smearing is for line collimation as in Kratky SAXS cameras and for point collimation for SANS. For SANS the resolution smearing a la Pedersen is realised (resFunct).
For desmearing the Lake algorithm is an iterative procedure to desmear smeared data. We follow here the improvements according to Vad using a convergence criterion and smoothing.
As references the waterXray scattering and a AgBe reference spectrum are available.
For form factors and structure factors see the respective modules.
5.1. Smear/desmear¶
smear (data, beamProfile, **kwargs) |
Smearing data for line-collimated SAXS (Kratky camera) or as point collimation SANS/SAXS. |
desmear (Ios, beamProfile[, maxIterations, …]) |
Desmearing according to Lake algorithm with a convergence criterion to stop recursion at best desmearing. |
resFunct (S[, collDist, collAperture, …]) |
Resolution smearing of small angle scattering for SANS or SAXS according to Pedersen for radial averaged data. |
prepareBeamProfile ([data]) |
Prepare beam profile parameter from Beam Profile measurement or according to given parameters. |
getBeamWidth (empty[, minmax, showfit]) |
Fit Gaussian to beam width of measurement. |
plotBeamProfile (beam[, p]) |
Plots beam profile and weight function according to parameters in beam. |
waterXrayScattering ([composition, T, units]) |
Absolute scattering of water with components (salt, buffer) at Q=0 as reference for X-ray. |
AgBeReference (q, wavelength[, n, ampn, …]) |
The scattering intensity expected from AgBe as a reference for wavelength calibration. |
5.2. Housekeeping¶
readpdh (pdhFileName) |
Opens and reads a SAXS data file in the .pdh (Primary Data Handling) format. |
autoscaleYinoverlapX (dataa[, key, keep]) |
Scales elements of data to have same mean .Y value in the overlap region of .X . |
removespikesminmaxmethod (dataa[, order, …]) |
Takes a dataset and removes single spikes from data by substitution with spline. |
removespikes (dataa[, xmin, xmax, medwindow, …]) |
Takes a dataset and removes single spikes. |
locatefiles (pattern[, root]) |
Locate all files matching supplied filename pattern in and below supplied root directory. |
copyfiles (pattern[, root, destination, link]) |
Copies all files matching pattern in tree below root to destination directory |
addXMLParameter (dataa) |
Adds the parameters stored in xml part of a .pdh file |
moveSAXSPACE (pattern[, root, destination, …]) |
Read SAXSPACE .pdh files and removes spikes by removespikes. |
This module allows smearing/desmearing of SAXS/SANS data.
Smearing is for line collimation as in Kratky SAXS cameras and for point collimation for SANS. For SANS the resolution smearing a la Pedersen is realised (resFunct).
For desmearing the Lake algorithm is an iterative procedure to desmear smeared data. We follow here the improvements according to Vad using a convergence criterion and smoothing.
As references the waterXray scattering and a AgBe reference spectrum are available.
For form factors and structure factors see the respective modules.
-
jscatter.sas.
AgBeReference
(q, wavelength, n=array([1, 2, 3, 4, 5, 6, 7, 8, 9]), ampn=[1, 1, 1, 1, 1, 1, 1, 1, 1, 1], domainsize=100, udw=0.1, asym=0, lg=1)[source]¶ The scattering intensity expected from AgBe as a reference for wavelength calibration.
The intensities assume a d-spacing of 5.8378 nm and a reduction of the intensity as q**-2. The domain size determines the width according to Scherrer equation [R188]. The first peak is at 1.076 1/nm. The result needs to be convoluted with the instrument resolution by resFunct or smear.
Parameters: q : array
wavevector
wavelength : float
wavelength
n : array of int
order of the peaks
ampn : list of float
amplitudes of the peaks
domainsize : float
Domainsize of AgBe crystals in nm. default 100 nm as is given in [R187].
udw : float
displacement u in Debye Waller factor exp(-u**2*q**2/3)
asym : float
factor asymmetry in Voigt function describing the peaks
lg : float
Lorenzian/gaussian fraction of both FWHM, describes the contributions of gaussian and lorenzian shape. See Voigt for details.
Returns: dataArray
References
[R187] (1, 2) T. C. Huang, H. Toraya, T. N. Blanton and Y. Wu X-ray Powder Diffraction Analysis of Silver Behenate, a Possible Low-Angle Diffraction Standard J. Appl.Cryst.(1993).26,180-184 [R188] (1, 2) Patterson, A. The Scherrer Formula for X-Ray Particle Size Determination Phys. Rev. 56 (10): 978–982 (1939) doi:10.1103/PhysRev.56.978.
-
jscatter.sas.
addXMLParameter
(dataa)[source]¶ Adds the parameters stored in xml part of a .pdh file
eg. in SAXSPACE .pdh files
Takes all lines with a first non whithespace ‘<’ as input from the original file
-
jscatter.sas.
autoscaleYinoverlapX
(dataa, key=None, keep='lowest')[source]¶ Scales elements of data to have same mean .Y value in the overlap region of .X .
Parameters: dataa : dataList
data to scale
key : string
Data are grouped into unique values of attribute key before scaling. E.g. to do it for a series of concentrations for each concentration individually.
keep : default ‘l’
If ‘l’ the lowest X are kept and higher X are scaled successively to next lower X. Anything else highest X are kept and other are scaled to next higher.
Returns: dataList
new scaled dataList
Notes
First data are sorted along .X.mean() scaling value is stored in .autoscalefactor
-
jscatter.sas.
copyfiles
(pattern, root='.', destination='copy', link=False)[source]¶ Copies all files matching pattern in tree below root to destination directory
Parameters: pattern : file pattern
Pattern used in fnmatch.filter
root : directory, default is os.curdir
Directory where to start
destination : dirname
Destination
link : bool
If True links are created.
-
jscatter.sas.
desmear
(Ios, beamProfile, maxIterations=15, windowsize=3, **kwargs)[source]¶ Desmearing according to Lake algorithm with a convergence criterion to stop recursion at best desmearing.
The convergence criterion reaches a minimum and the Lake algorithm [R189] is stopped as described in [R190]. In each step a smoothing based on the ratio desmeared/observed as described in [R190] is used (point average with windowsize 3).
Parameters: Ios : dataArray
Original smeared data
beamProfile : dataArray
Beam profile as prepeared with prepareBeamProfile
maxIterations : int, default=15
maximum iterations to stop
windowsize : int , default=3
Window size for smoothing in each step of desmearing (running average).
Returns: dataArray
References
[R189] (1, 2) Lake, J. A. (1967). Acta Cryst. 23, 191–194. [R190] (1, 2, 3) Comparison of iterative desmearing procedures for one-dimensional small-angle scattering data Vad and Sager, J. Appl. Cryst. (2011) 44,32-42
-
jscatter.sas.
fitdarkcurrent
(darkfile)[source]¶ Fits dark current with 5th order polynom + cosine.
This is dangerous as the darkcurrent has a noise that is not random but depends on count time.
-
jscatter.sas.
getBeamWidth
(empty, minmax=[-0.03, 0.03], showfit=False)[source]¶ Fit Gaussian to beam width of measurement.
Parameters: empty : dataArray
Measurement with the transmitted beam included. This is used within min max to fit a Gaussian.
minmax : [float,float]
Minimum and maximum for fitting
showfit : bool
Show the fit result
Returns: float as the beam width
-
jscatter.sas.
locatefiles
(pattern, root='.')[source]¶ Locate all files matching supplied filename pattern in and below supplied root directory.
Parameters: pattern : file pattern
Pattern used in fnmatch.filter
root : directory, default is os.curdir
Directory where to start
Returns: file list
-
jscatter.sas.
moveSAXSPACE
(pattern, root='./', destination='./despiked', skip='BeamProfile', medwindow=5, SGwindow=5, sigma=0.2, order=2)[source]¶ Read SAXSPACE .pdh files and removes spikes by removespikes.
This is mainly for use at JCNS SAXSPACE with CCD camera as detector :-))))
Parameters: pattern : string
Search pattern for filenames
destination : string
Where to save the files
skip : string
Words in filename to skip the file
medwindow : odd integer
Window size of scipy.signal.medfilt
SGwindow : odd int, None
Savitzky-Golay filter see scipy.signal.savgol_filter
order : int
Polynominal order of scipy.signal.savgol_filter
sigma : float
Deviation factor of eY If datapoint-median> sigma*std its a spike
Notes
Default values are adjusted to typical SAXSPACE measurement SAXS not to reduce peaks and transmission peak and overshot region excluded sigma =0.2 because the measurement error seems to be to big
-
jscatter.sas.
plotBeamProfile
(beam, p=None)[source]¶ Plots beam profile and weight function according to parameters in beam.
Parameters: beam
beam with parameters
-
jscatter.sas.
prepareBeamProfile
(data=None, **kwargs)[source]¶ Prepare beam profile parameter from Beam Profile measurement or according to given parameters.
Parameters: data : dataArray,’trapez’,’SANS’
- Contains measured beam profile, explicit Gaussian width list or type ‘SANS’, ‘trapz’.
- dataArray Line collimation as measured can be given and will be smoothed and made symmetric.
- dataArray with explicit given Gaussian width missing .X values will be interpolated.
- ‘trapez’ : Line collimation with trapezoidal parameters a, b, bxw, dIW.
- ‘SANS’ : Smearing a la Pedersen; see resFunct for parameters
collDist,collAperture,detDist,sampleAperture,wavelength,wavespread,dpixelWidth,dringwidth,zeroextrapolfunc :
Parameters as described in resFunct for SANS These are determined from the experimental setup.
a : float
Larger full length of trapezoidal profile in detector q units
b : float
Shorter full length of trapezoidal profile in detector q units If a=b ==> a=a*(1+1e-7), b=b*(1-1e-7)
bxw : float
Beam half-width (at half-height) that defines the Gaussian to describe the width profile in detector q units. See getBeamWidth to fit to a measurement.
dIW : float
Detector slit width in detector q units. Length on detector to integrate parallel to beam length for line collimation.
wavelength : float,
Wavelength in nm default 0.155418 nm for SAXS 0.5 for SANS
detDist : float, default 305.3558
Detectordistance in units mm Default is SAXS detector distance of SAXSpace
explicit : int
For explicit given Gaussian width the index of the column with the width. For merged dataFiles of KWS2@MLZ this is the forth column with index 3.
Returns: beam profile as dataArray
Notes
- For measured beam profiles parameters a,b are determined from the flanks for trapezoidal profile.
- Detector q units are equivalent to the pixel distance as expressed in a corrected measurement.
- For ‘explicit’ Gaussian width a SANS measurement as on KWS2 can be used which has sigma as 4th column. Missing values are interpolated.
Examples
# use as # prepare measured beamprofile mbeam = js.sas.prepareBeamProfile(beam, bxw=0.01, dIW=1.) # prepare profile with trapezoidal shape (a,b are fitted above) tbeam = js.sas.prepareBeamProfile('trapz', a=mbeam.a, b=mbeam.b, bxw=0.01, dIW=1) # prepare profile SANS (missing parameters get defaults) Sbeam = js.sas.prepareBeamProfile('SANS', detDist=2000,wavelength=0.4,wavespread=0.15) # prepare profile with explicit given Gaussian width in column 3 as e.g. KWS2@MLZ Gbeam = js.sas.prepareBeamProfile(measurement,explicit=3)
-
jscatter.sas.
readpdh
(pdhFileName)[source]¶ Opens and reads a SAXS data file in the .pdh (Primary Data Handling) format.
If data contain X values <0 it is assumed that the primary beam is included as for SAXSpace instruments. In this case the primary beam with max value as .transmission and the peak center as centerTransmissionPeak are extracted.
Parameters: pdhFileName : string
file name
Returns: dataArray
Notes
PDH format used in the PCG SAXS software suite developed by the Glatter group at the University of Graz, Austria. This format is described in the appendix of the PCG manual (page 123 of the 2005 version). In the PDH format, lines 1-5 contain header information, followed by the SAXS data. A header looks like this:
saxs082511water SAXS 470 1 0 0 0 0 0 0 0.000000E+00 0.264500E+03 0.000000E+00 0.100000E+01 0.154200E+00 0.818682E+02 0.306286E-01 0.181261E-02 0.154053E+00 0.100000E+01 .....
- line 3:
- 0: number of points
- line 4:
- 0:
- 1: detector distance => SAXSpace
- 3: normalization factor; not used
- 4: wavelength => SAXSpace
- line 5:
- 0:
- 1: beam-length profile parameter a
- 2: beam-length profile parameter b
- 3: half-width of the beam-width profile (at half-height)
- 4: detector slit length (equivalent to width of integration area)
- 5: q units scale: 1.0 for A-1, 10.0 for nm-1
Additional xml parameter as in the SAXSPACE format can be extracted to attributes by addXMLParameter.
-
jscatter.sas.
removespikes
(dataa, xmin=None, xmax=None, medwindow=5, SGwindow=None, sigma=0.2, SGorder=2)[source]¶ Takes a dataset and removes single spikes.
A median filter is used to find single spikes. If abs(data.Y-medianY)/data.eY>sigma then the medianY value is used. If SGwindow!=None Savitzky-Golay filtered values are used. If sigma is 0 then new values (median or Savitzky-Golay filtered) are used everywhere.
Parameters: dataa : dataArray
dataset with eY data
medwindow : odd integer
window size of scipy.signal.medfilt
SGwindow : odd int, None
Savitzky-Golay filter see scipy.signal.savgol_filter without the spikes; window should be smaller than instrument resolution
order : int
polynominal order of scipy.signal.savgol_filter needs to be smaller than SGwindow
SGsigma : float
relative deviation from eY if datapoint-median> sigma*eY its a spike
Returns: Filtered and smoothed dataArray
-
jscatter.sas.
removespikesminmaxmethod
(dataa, order=7, sigma=2, nrepeat=1, removePoints=None)[source]¶ Takes a dataset and removes single spikes from data by substitution with spline.
Find minima and maxima of data including double point spikes; no 3 point spikes scipy.signal.argrelextrema with np.greater and np.less are used to find extrema
Parameters: dataa : dataArray
Dataset with eY data.
order : int
Number of points see scipy.signal.argrelextrema. Distance between extrema.
sigma : float
Deviation factor from std dev; from eY. If datapoint -spline> sigma*std its a spike.
nrepeat : int
Repeat the procedure nrepeat times.
removePoints : list of integer
Instrument related points to remove because of dead pixels. ‘JCNS’ results in a list for SAXSPACE at JCNS Jülich.
Returns: data with spikes removed
-
jscatter.sas.
resFunct
(S, collDist=8000.0, collAperture=10, detDist=8000.0, sampleAperture=10, wavelength=0.5, wavespread=0.2, dpixelWidth=10, dringwidth=1, zeroextrapolfunc=2, highQextrapolfunc=2)[source]¶ Resolution smearing of small angle scattering for SANS or SAXS according to Pedersen for radial averaged data.
I(q0)= Integral{(R(q,q0)*S(q)}dq with Kernel R(q,q0) of equ. 33 in [R191] including wavelength spread, finite collimation and detector resolution. Default parameters are typical for a SANS machine like KWS2@JCNS with rectangular apertures. Low Q can be extrapolated as power law or Guinier like or constant.
Parameters: S : array like
dataArray with X. and .Y theoretical Scattering function q in nm^-1 .Y can be arbitrary unit
collDist : float, default 8000
collimation distance in mm
collAperture : float, default 10
collimation rectangular aperture size in mm
detDist : float, default 8000
detector distance in mm
sampleAperture : float, default 10
sample rectangular aperture size in mm
wavelength : float, default 0.5
wavelength in nm
wavespread : float, default 0.1
FWHM wavelengthspread dlambda/lambda
dpixelWidth : float, default 10
Detector pixel width in mm
dringwidth : integer, default 1
number of pixel for averaging
zeroextrapolfunc : , default 2
- Type of extrapolation at low X edge for better handling of the border:
- guinier : Low X data are log scaled, then X**2 extrapolated as Guinier like extrapolation.
- float : Power law extrapolation of low X e.g. -4 for X**-4 for Porod scaling.
- None : A constant value as Y(X.min()).
highQextrapolfunc : float, default 2
- Type of extrapolation at high X edge for better handling of the border:
- float : Power law extrapolation of low X e.g. -4 for X**-4 for Porod scaling.
- None : A constant value as Y(X.max()).
Returns
——-
dataArray
columns [‘wavevector; smeared scattering; unsmeared scattering; half width smearing function’]
Notes
- HalfWidthSmearingFunction is the FWHM the Gaussian used for smearing including all effects.
- The resolution is assumed to be equal in direction parallel and perpendicular to q on a 2D detector as described in chap. 2.5 in [R191].
- We neglect additional smearing due to radial averaging (last paragraph in chap 2.5 of [R191]).
- Defaults correspond to a typical medium resolution measurement.
- zeroextrapolfunc, highQextrapolfuncallow extrapolation at both edges to reduce edge effects. The best values depends on the measured signal shape at the edge and may change. The optimal way is to calculate the used model for the whole Q range, smear it and prune to the needed range. This is demonstated in example 2.
- An example for SANS fitting with resFunc is given in example_SANSsmearing.py.
References
[R191] (1, 2, 3, 4, 5) Analytical Treatment of the Resolution Function for Small-Angle Scattering JAN SKOV PEDERSEN, DORTHE POSSELTAND KELL MORTENSEN J. Appl. Cryst. (1990). 23, 321-333 Examples
Reproducing Table 1 of [R191]
import jscatter as js q=js.loglist(0.1,10,500) S=js.ff.sphere(q,6) # this is the direct call of resFunc, use smear instead as shown in next example Sr=js.sas.resFunct(S, collDist=2000.0, collAperture=20, detDist=2000.0, sampleAperture=10, wavelength=0.5, wavespread=0.2,dpixelWidth=0,dringwidth=0) # plot it p=js.grace() p.plot( S,sy=[1,0.3],li=1,legend='sphere') p.plot( Sr,sy=[2,0.3,2],li=2,legend='smeared sphere') p.plot(Sr.X,Sr[-1],li=4,sy=0,legend='FWHM in nm\S-1 ') p.yaxis(min=1e-3,scale='l',charsize=1.5,label='I(q) / a.u.',tick=[10,9]) p.yaxis(min=1e-1,tick=[10,9]) p.xaxis(scale='l',charsize=1.5,label='q / nm\S-1',tick=[10,9]) p.legend(x=0.8,y=5e5)
Example 2:
# smear model over full range and interpolate to needed data # this is the best way to smear a model for fitting, but is not possible for desmearing meas=js.dA('measureddata.dat') # load data # define profile resol2m = js.sas.prepareBeamProfile('SANS', detDist=2000,collDist=2000.,wavelength=0.4,wavespread=0.15) q=np.r_[0.01:5:0.01] # or np.r_[0:meas.X.min():0.01,meas.X,meas.X.max():meas.X.max()*2:0.1] # calc model temp=js.ff.ellipsoid(q,2,3) # smear it smearedmodel=js.sas.smear(temp,resol2m).interpolate(X=meas.X)
-
jscatter.sas.
resFunctExplicit
(S, beamprofile, zeroextrapolfunc=2, highQextrapolfunc=2)[source]¶ Resolution smearing of small angle scattering for SANS or SAXS according to explict given Gaussian width.
I(q0)= Integral{(R(q,q0)*S(q)}dq with Gaussian kernel R(q,q0). Low Q can be extrapolated as power law or Guinier like or constant.
Parameters: S : array like
dataArray with X. and .Y theoretical Scattering function q in nm^-1 .Y can be arbitrary unit
beamProfile : beamProfile ‘explicit’
Beam profile as prepared from prepareBeamProfile ‘explicit’
zeroextrapolfunc : , default 2
- Type of extrapolation at low X edge for better handling of the border:
- guinier : Low X data are log scaled, then X**2 extrapolated as Guinier like extrapolation.
- float : Power law extrapolation of low X e.g. -4 for X**-4 for Porod scaling.
- None : A constant value as Y(X.min()).
highQextrapolfunc : float, default 2
- Type of extrapolation at high X edge for better handling of the border:
- float : Power law extrapolation of low X e.g. -4 for X**-4 for Porod scaling.
- None : A constant value as Y(X.max()).
Returns
——-
dataArray
columns [‘wavevector; smeared scattering; unsmeared scattering; half width smearing function’]
Notes
- HalfWidthSmearingFunction is the FWHM the Gaussian used for smearing including all effects.
-
jscatter.sas.
smear
(data, beamProfile, **kwargs)[source]¶ Smearing data for line-collimated SAXS (Kratky camera) or as point collimation SANS/SAXS.
The full resolution for SANS is described in resFunct.
Parameters: data : dataArray
Data to be smeared
beamProfile : beamProfile or ‘trap’, ‘SANS’, ‘explicit’, dataArray
Beam profile as prepared from prepareBeamProfile or type as ‘trapezoidal’, ‘SANS’,’explicit’ or a measured beam profile as dataArray for line collimation. Measured Profile is treated by prepareBeamProfile.
kwargs :
See prepareBeamProfile for kwargs.
Returns: dataArray
Notes
- If data has attributes a, b, dIW, bxw, detDist these are used, if not given in function call.
- If wavelength is missing in data a default of 0.155418 nm for Xray K_alpha line is assumed. For SANS 0.6 A.
Examples
# use as # prepare measured line collimation beamprofile mbeam = js.sas.prepareBeamProfile(beam, bxw=0.01, dIW=1.) # prepare profile with trapezoidal shape (a,b can be fitted above) tbeam = js.sas.prepareBeamProfile('trapz', a=mbeam.a, b=mbeam.b, bxw=0.01, dIW=1) # prepare profile SANS (missing parameters get defaults, see resFunct) Sbeam = js.sas.prepareBeamProfile('SANS', detDist=2000,wavelength=0.4,wavespread=0.15) # prepare profile with explicit given Gaussian width in column 3 as e.g. KWS2@JCNS Gbeam = js.sas.prepareBeamProfile(measurement,explicit=3) # smear datasm= js.sas.smear(data,mbeam) datast= js.sas.smear(data,tbeam) datasS= js.sas.smear(data,Sbeam) datasG= js.sas.smear(data,Gbeam)
-
jscatter.sas.
waterXrayScattering
(composition='h2o1', T=293, units='mol')[source]¶ Absolute scattering of water with components (salt, buffer) at Q=0 as reference for X-ray.
According to [R193] a buffer of water with components might be used. Ions need to be given separatly as [‘55.51h2o1’,‘0.15Na’,‘0.15Cl’] for 0.15 M NaCl solution. It is accounted for the temperature dependence of water density and compressibility.
Parameters: composition : string
Buffer composition as in scatteringLengthDensityCalc give dissosiated ions separatly as [‘1Na’,‘1Cl’] with concentration in mol prepended the additional scattering as ionic liquid of the ions in water is taken into account see [R193] mass in g; 1000g water are 55.508 mol
T : float
temperature in °K
units : ‘mol’
anything except ‘mol’ prepended unit is mass fraction ‘mol’ prepended units is mol and mass fraction is calculated as
e.g. 1l Water with 123mmol NaCl [‘55.508H2O1’,‘0.123Na1Cl1’]
Returns: float absolute scattering length in Units 1/cm
Notes
with
water density
compressibility
number of electrons per water molecule
cross section of electron in nm
Boltzman constant
concentration component i
number of electrons per molecule component i in Mol
is done for all ions separately if given
References
[R192] SAXS experiments on absolute scale with Kratky systems using water as a secondary standard Doris Orthaber et al. J. Appl. Cryst. (2000). 33, 218±225 [R193] (1, 2, 3) A high sensitivity pinhole camera for soft condensed matter T. Zemb, O. Tache, F. Né, and O. Spalla, J. Appl. Crystallogr. 36, 800 (2003).