5. smallanglescattering (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[, NIterations, …]) |
Desmearing according to Lake algorithm with posibility 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 from Beam Profile measurement or according to given parameters. |
getBeamWidth (empty[, minmax, show]) |
Extract primary beam of empty cell or buffer measurement. |
plotBeamProfile (beam[, p]) |
Plots beam profile and weight function according to parameters in beam. |
transmissionCorrection (data, dark[, …]) |
Subtract dark current, find primary beam, get transmission and normalize by transmission and exposure time. |
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 (data) |
Adds the parameters stored in xml part of a .pdh file as eg. |
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.smallanglescattering.
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 [R196]. 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 [R195].
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
[R195] (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 [R196] (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.smallanglescattering.
addXMLParameter
(data)[source]¶ Adds the parameters stored in xml part of a .pdh file as eg. in SAXSPACE .pdh files.
Parameters: data : dataArray
Already read pdh file. XML content is found in comments of the read files and starts with ‘<’.
-
jscatter.smallanglescattering.
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.smallanglescattering.
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.smallanglescattering.
desmear
(Ios, beamProfile, NIterations=-15, windowsize=4, qmax=4, **kwargs)[source]¶ Desmearing according to Lake algorithm with posibility to stop recursion at best desmearing.
For negative NIterations the iterations are stopped if a convergence criterion reaches a minimum as described by Vad [R198]. In each step a smoothing based on the ratio desmeared/observed as described in [R198] is used (point average with windowsize).
Parameters: Ios : dataArray
Original smeared data
beamProfile : dataArray
Beam profile as prepeared with prepareBeamProfile
NIterations : int, default=-15
Number of iterations to stop. Negative values indicate to use the convergence criterion described by Vad [R198] and abs(NIterations) gives the maximum iterations to stop.
qmax : float, default=4
Maximum in scattering vector q up to where the convergence criterion is evaluated. This reduces the influence of the noise at larger a.
windowsize : odd int , default=4
Window size for smoothing in each step of desmearing (running average).
Returns: dataArray
References
[R197] Lake, J. A. (1967). Acta Cryst. 23, 191–194. [R198] (1, 2, 3, 4) Comparison of iterative desmearing procedures for one-dimensional small-angle scattering data Vad and Sager, J. Appl. Cryst. (2011) 44,32-42
-
jscatter.smallanglescattering.
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 used detector and count time.
-
jscatter.smallanglescattering.
getBeamWidth
(empty, minmax='auto', show=False)[source]¶ Extract primary beam of empty cell or buffer measurement.
The primary beam is searched and cut between the next minima found, then normalized. Additionally a Gaussian fit is done and hwhm is included in result profile.
Parameters: empty : dataArray
Empty cell measurement with the transmitted beam included.
minmax : ‘auto’,[float,float]
Automatic or interval for search of primary beam. E.g. [-0.03,0.03] allow for explicitly setting the interval.
show : bool
Show the fit result
Returns: dataArray with beam width profile
.sigma sigma of fit with Gaussian .hwhm half width half maximum
-
jscatter.smallanglescattering.
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.smallanglescattering.
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.
-
jscatter.smallanglescattering.
plotBeamProfile
(beam, p=None)[source]¶ Plots beam profile and weight function according to parameters in beam.
Parameters: beam
beam with parameters
-
jscatter.smallanglescattering.
prepareBeamProfile
(data=None, **kwargs)[source]¶ Prepare beam profile 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 for each Q values, missing 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,dataArray
Beam width profile. Use getBeamWidth to cut the primary beam and fit a Gaussian. A float describes the beam half-width at half maximum (hwhm of Gaussian). If bxw is the profile prepared by getBeamWidth the experiemntal profile is used.
dIW : float
Detector slit width in detector q units. Length on detector to integrate parallel to beam length for line collimation. On my SAXspace this is 1.332 as given in the header if the file.
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.smallanglescattering.
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
Alternativly the files can be read ignoring the information in the header (it is stored in the comments)
data=dA(pdhFileName,lines2parameter=[2,3,4])
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 (below from version 4.05.12 page 159).
In the PDH format, lines 1-5 contain header information, followed by the SAXS data.
line 1: format A80 -> description line 2: format 16(A4,1X) -> description in 16x4 character groups (1X = space separated) line 3: format 8(I9,1X) -> 8 integers (1X = space separated) line 4: format 5(E14.6,1X) -> 5 float (1X = space separated) line 5: format 5(E14.6,1X -> 5 float(1X = space separated) line 6+ format 3(E14.6,1X) - SAXS data x,y,error (1X = space separated)
- with:
- line 3 field 0 : number of points
- line 4 field 4 : normalization constant (default 1, never zero!!)
Anything else can have different meanings.
- The SAXSpace and SAXSess (AntonPaar) format add:
- line 4 field 2 : detector distance in mm
- line 4 field 5 : wavelength
- line 5 field 2 : detector slit length (equivalent to width of integration area) in q units
Additional xml parameter as in the SAXSPACE format appended can be extracted to attributes by addXMLParameter. Mainly this is “Exposure” time.
Example data for SAXSpace
<emptyline> SAXS BOX 2057 0 0 0 0 0 0 0 0.000000E+00 3.052516E+02 0.000000E+00 1.000000E+00 1.541800E-01 0.000000E+00 1.332843E+00 0.000000E+00 0.000000E+00 0.000000E+00 -1.033335E-01 2.241656E+03 1.024389E+00 -1.001430E-01 2.199972E+03 1.052537E+00 ...
-
jscatter.smallanglescattering.
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.smallanglescattering.
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.smallanglescattering.
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 [R199] 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 [R199].
- We neglect additional smearing due to radial averaging (last paragraph in chap 2.5 of [R199]).
- 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
[R199] (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 [R199]
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.smallanglescattering.
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). E.g. for merged dataFiles of KWS2@MLZ the explicit width is given in the 4th column.
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.smallanglescattering.
smear
(data, beamProfile, **kwargs)[source]¶ Smearing data for line-collimated SAXS (Kratky camera) or as point collimation SANS/SAXS.
The full resolution for point collimation SAXS/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.
- During smearing for Kratky camera an integration over the beam width and beam length are performed. In this integration q_{w,l}= ((q+q_{w})^2)+q_{l}^2)^{1/2} is used with q_{w} along the beam width and q_{l} along the beam length. in regions q_{w,l} > max(q_{data}) we estimate the measured scattering intensity by the mean of the last 10 points of the measured spectra to allow for a maximum in q range. The strictly valid q range can be estimated by calculating q_{x,y} < max(q) with 2 times the used beam width and beam length. As the smearing for larger q has no real effect the estimate might be still ok.
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.smallanglescattering.
transmissionCorrection
(data, dark, emptybeam=None, edge=0.03, exposure=None)[source]¶ Subtract dark current, find primary beam, get transmission and normalize by transmission and exposure time.
For measurements including the primary beam from a semitransparent beamstop as from SAXSpace. Transmission is the maximum of the primary beam peak after dark subtraction. Allows easier comaprison of measurements with aged source (primary intensity change).
Parameters: data : dataArray
A measurement from a SAXSpace instrument read by js.dA(‘filename.pdh’,lines2parameter=[2,3,4])
dark : dataArray
Dark current measurement.
emptybeam : dataArray,float, default=None
- Empty beam measurement or count rate at primary peak maximum of emty peak.
- dark will be subtracted.
- If not provided the empty beam transmission is set to 1.
- If float the empty beam is not subtracted and only transmission is calculated as T=\frac{I(q=0)_{sample}}{I(q=0)_{emptybeam}}.
edge : float, default 0.03
Wavevector value below beam stop edge. The primary beam is searched below this value.
exposure : float, default None
Esposure time in unit ‘s’. If not given the xml description at the end of the file is examined.
Returns: dataArray with corrected data as \frac{I_S-I_{dark}}{T_S} (see Notes)
Notes
Files from SAXSpace instrument (.pdh) can be read by
js.dA('filename.pdh',lines2parameter=[2,3,4])
Correction
Brulet at al [R200] describe the data correction for SANS, which is in principle also valid for SAXS, if incoherent contributions are neglected.
The difference is, that SAXS has typical transmission around ~0.3 for 1mm water sample in quartz cell due to absorption, while in SANS typical values are around ~0.9 for D2O. Larger volume fractions in the sample play a more important rule for SANS as hydrogenated ingredients reduce the transmission significantly, while in SAXS still the water and the cell (quartz) dominate.
One finds for a sample inside of a container with thicknesses (z) for sample, buffer (solvent), empty cell and emty beam measurement (omitting the overall q dependence):
I_s = \frac{1}{z_S}\big((\frac{I_S-I_{dark}}{T_S}-I_{b}T_S\big) -\big(\frac{I_{EC}-I_{dark}}{T_{EC}}-I_{b}T_{EC})\big) - \frac{1}{z_B}\big((\frac{I_B-I_{dark}}{T_B}-I_{b}T_B\big) -\big(\frac{I_{EC}-I_{dark}}{T_{EC}}-I_{b}T_{EC})\big)
- where
- I_s is the interesting species
- I_S is the sample of species in solvent (buffer)
- I_B is the pure solvent (describing a constant background)
- I_{dark} is the dark current measurement
- I_b is the empty beam measurement
- I_{EC} is the empty cell measurement
- z_x corresponding sample thickness
- T_x corresponding transmission
The reccuring pattern \big((\frac{I-I_{dark}}{T}-I_{b}T\big) shows that the the beam tail (border of primary beam not absorbed by the beam stop) is attenuated by the corresponding sample.
For equal sample thickness z the empty beam is included in subtraction of I_B:
I_s = \frac{1}{z} \big((\frac{I_S-I_{dark}}{T_S}-I_{b}T_S) - (\frac{I_B-I_{dark}}{T_B}-I_{b}T_B)\big)
The simple case
If the transmissions are nearly equal as for e.g. protein samples with low concentration (T_S \approx T_B) we only need to subtract the transmission and dark current corrected buffer measurement from the sample.
I_s = \frac{1}{z} \big((\frac{I_S-I_{dark}}{T_S}) - (\frac{I_B-I_{dark}}{T_B}\big)
Higher accuracy for large volume fractions
For larger volume fractions \Phi the transmission might be different and we have to take into account that only 1-\Phi of solvent contributes to I_S. We may incorporate this in the sense of an optical density changing the effective thickness \frac{1}{z_B}\rightarrow\frac{1-\Phi}{z_B} resulting in different thicknesses z_S \neq z_B
Transmission
The transmission is measured as the ratio T=\frac{I(q=0)_{sample}}{I(q=0)_{emptybeam}} with I(q=0) as the primary beam intensity.
If the primary beam tail is neglected in the above equation I(q=0)_{emptybeam} only gives a common scaling factor and can be omitted if arbitrary units are used. Alternatively one can scale to the EC transmission with T_{EC}=1 For absolute calibration the same needs to be done. One finds T_{sample in cell}=T_{empty cell}T_{sample without cell}.
References
[R200] (1, 2) Improvement of data treatment in small-angle neutron scattering Brûlet et al Journal of Applied Crystallography 40, 165-177 (2007)
-
jscatter.smallanglescattering.
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 [R202] 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 [R202] 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 mass=[mol] mass_{molecule} e.g. 1l Water with 123mmol NaCl [‘55.508H2O1’,‘0.123Na1Cl1’]
Returns: float absolute scattering length in Units 1/cm
Notes
I(0)=(\sigma_{water}^2f_e^2 n_{ew}^2 k_B T \chi + \sum_{ci} n_i N_A 1000 n_{ei}^2 f_e^2 )/100
with
\sigma_{water} water density \chi compressibility n_{ew} number of electrons per water molecule f_e cross section of electron in nm k_B Boltzman constant n_i concentration component i n_{ei} number of electrons per molecule component i in Mol \sum_{ci} is done for all ions separately if givenReferences
[R201] 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 [R202] (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).