xrayutilities.analysis package¶
Submodules¶
xrayutilities.analysis.line_cuts module¶
- xrayutilities.analysis.line_cuts.get_arbitrary_line(qpos, intensity, point, vec, npoints, intrange)[source]¶
extracts a line scan from reciprocal space map data along an arbitrary line defined by the point ‘point’ and propergation vector ‘vec’. Integration of the data is performed in a cylindrical volume along the line. This function works for 2D and 3D input data!
- Parameters:
- qposlist of array-like objects
arrays of x, y (list with two components) or x, y, z (list with three components) momentum transfers
- intensityarray-like
2D or 3D array of reciprocal space intensity with shape equal to the qpos entries
- pointtuple, list or array-like
point on the extraction line (2 or 3 coordinates)
- vectuple, list or array-like
propergation vector of the extraction line (2 or 3 coordinates)
- npointsint
number of points in the output data
- intrangefloat
radius of the cylindrical integration volume around the extraction line
- Returns:
- qpos, qintndarray
line scan coordinates and intensities
- used_maskndarray
mask of used data, shape is the same as the input intensity: True for points which contributed, False for all others
Examples
>>> qcut, qint, mask = get_arbitrary_line([qx, qy, qz], inten, ... (1.1, 2.2, 0.0), ... (1, 1, 1), 200, 0.1)
- xrayutilities.analysis.line_cuts.get_omega_scan(qpos, intensity, cutpos, npoints, intrange, **kwargs)[source]¶
extracts an omega scan from reciprocal space map data with integration along either the 2theta, or radial (omega-2theta) direction. The coplanar diffraction geometry with qy and qz as diffraction plane is assumed. This is consistent with the coplanar geometry implemented in the HXRD-experiment class.
This function works for 2D and 3D input data in the same way!
- Parameters:
- qposlist of array-like objects
arrays of y, z (list with two components) or x, y, z (list with three components) momentum transfers
- intensityarray-like
2D or 3D array of reciprocal space intensity with shape equal to the qpos entries
- cutpostuple or list
y/z-position or x/y/z-position at which the line scan should be extracted. this must be have two entries for 2D data (z-position) and a three for 3D data
- npointsint
number of points in the output data
- intrangefloat
integration range in along intdir in degree. data will be integrated from -intrange .. +intrange
- intdir{‘2theta’, ‘radial’}, optional
integration direction: ‘2theta’: scattering angle (default), or ‘radial’: omega-2theta direction.
- wlfloat or str, optional
wavelength used to determine angular integration positions
- Note:
Although applicable for any set of data, the extraction only makes sense when the data are aligned into the y/z-plane.
- Returns:
- om, omintndarray
omega scan coordinates and intensities
- used_maskndarray
mask of used data, shape is the same as the input intensity: True for points which contributed, False for all others
Examples
>>> omcut, omcut_int, mask = get_omega_scan([qy, qz], inten, [2.0, 5.0], ... 250, intrange=0.1)
- xrayutilities.analysis.line_cuts.get_qx_scan(qpos, intensity, cutpos, npoints, intrange, **kwargs)[source]¶
extracts a qx scan from 3D reciprocal space map data with integration along either, the perpendicular plane in q-space, omega (sample rocking angle) or 2theta direction. For the integration in angular space (omega, or 2theta) the coplanar diffraction geometry with qy and qz as diffraction plane is assumed. This is consistent with the coplanar geometry implemented in the HXRD-experiment class.
- Parameters:
- qposlist of array-like objects
arrays of x, y, z (list with three components) momentum transfers
- intensityarray-like
3D array of reciprocal space intensity with shape equal to the qpos entries
- cutpostuple/list
y/z-position at which the line scan should be extracted. this must be and a tuple/list with the qy, qz cut position
- npointsint
number of points in the output data
- intrangefloat
integration range in along intdir, either in 1/AA (q) or degree (‘omega’, or ‘2theta’). data will be integrated from -intrange .. +intrange
- intdir{‘q’, ‘omega’, ‘2theta’}, optional
integration direction: ‘q’: perpendicular Q-plane (default), ‘omega’: sample rocking angle, or ‘2theta’: scattering angle.
- wlfloat or str, optional
wavelength used to determine angular integration positions
- Note:
The angular integration directions although applicable for any set of data only makes sense when the data are aligned into the y/z-plane.
- Returns:
- qx, qxintndarray
qx scan coordinates and intensities
- used_maskndarray
mask of used data, shape is the same as the input intensity: True for points which contributed, False for all others
Examples
>>> qxcut, qxcut_int, mask = get_qx_scan([qx, qy, qz], inten, [0, 2.0], \ ... 250, intrange=0.01)
- xrayutilities.analysis.line_cuts.get_qy_scan(qpos, intensity, cutpos, npoints, intrange, **kwargs)[source]¶
extracts a qy scan from reciprocal space map data with integration along either, the perpendicular plane in q-space, omega (sample rocking angle) or 2theta direction. For the integration in angular space (omega, or 2theta) the coplanar diffraction geometry with qy and qz as diffraction plane is assumed. This is consistent with the coplanar geometry implemented in the HXRD-experiment class.
This function works for 2D and 3D input data in the same way!
- Parameters:
- qposlist of array-like objects
arrays of y, z (list with two components) or x, y, z (list with three components) momentum transfers
- intensityarray-like
2D or 3D array of reciprocal space intensity with shape equal to the qpos entries
- cutposfloat or tuple/list
x/z-position at which the line scan should be extracted. this must be a float for 2D data (z-position) and a tuple with two values for 3D data
- npointsint
number of points in the output data
- intrangefloat
integration range in along intdir, either in 1/AA (q) or degree (‘omega’, or ‘2theta’). data will be integrated from -intrange .. +intrange
- intdir{‘q’, ‘omega’, ‘2theta’}, optional
integration direction: ‘q’: perpendicular Q-plane (default), ‘omega’: sample rocking angle, or ‘2theta’: scattering angle.
- wlfloat or str, optional
wavelength used to determine angular integration positions
- Note:
For 3D data the angular integration directions although applicable for any set of data only makes sense when the data are aligned into the y/z-plane.
- Returns:
- qy, qyintndarray
qy scan coordinates and intensities
- used_maskndarray
mask of used data, shape is the same as the input intensity: True for points which contributed, False for all others
Examples
>>> qycut, qycut_int, mask = get_qy_scan([qy, qz], inten, 5.0, 250, \ ... intrange=0.02, intdir='2theta')
- xrayutilities.analysis.line_cuts.get_qz_scan(qpos, intensity, cutpos, npoints, intrange, **kwargs)[source]¶
extracts a qz scan from reciprocal space map data with integration along either, the perpendicular plane in q-space, omega (sample rocking angle) or 2theta direction. For the integration in angular space (omega, or 2theta) the coplanar diffraction geometry with qy and qz as diffraction plane is assumed. This is consistent with the coplanar geometry implemented in the HXRD-experiment class.
This function works for 2D and 3D input data in the same way!
- Parameters:
- qposlist of array-like objects
arrays of y, z (list with two components) or x, y, z (list with three components) momentum transfers
- intensityarray-like
2D or 3D array of reciprocal space intensity with shape equal to the qpos entries
- cutposfloat or tuple/list
x/y-position at which the line scan should be extracted. this must be a float for 2D data and a tuple with two values for 3D data
- npointsint
number of points in the output data
- intrangefloat
integration range in along intdir, either in 1/AA (q) or degree (‘omega’, or ‘2theta’). data will be integrated from -intrange/2 .. +intrange/2
- intdir{‘q’, ‘omega’, ‘2theta’}, optional
integration direction: ‘q’: perpendicular Q-plane (default), ‘omega’: sample rocking angle, or ‘2theta’: scattering angle.
- wlfloat or str, optional
wavelength used to determine angular integration positions
- Note:
For 3D data the angular integration directions although applicable for any set of data only makes sense when the data are aligned into the y/z-plane.
- Returns:
- qz, qzintndarray
qz scan coordinates and intensities
- used_maskndarray
mask of used data, shape is the same as the input intensity: True for points which contributed, False for all others
Examples
>>> qzcut, qzcut_int, mask = get_qz_scan([qy, qz], inten, 3.0, 200,\ ... intrange=0.3)
- xrayutilities.analysis.line_cuts.get_radial_scan(qpos, intensity, cutpos, npoints, intrange, **kwargs)[source]¶
extracts a radial scan from reciprocal space map data with integration along either the omega or 2theta direction. The coplanar diffraction geometry with qy and qz as diffraction plane is assumed. This is consistent with the coplanar geometry implemented in the HXRD-experiment class.
This function works for 2D and 3D input data in the same way!
- Parameters:
- qposlist of array-like objects
arrays of y, z (list with two components) or x, y, z (list with three components) momentum transfers
- intensityarray-like
2D or 3D array of reciprocal space intensity with shape equal to the qpos entries
- cutpostuple or list
y/z-position or x/y/z-position at which the line scan should be extracted. this must be have two entries for 2D data (z-position) and a three for 3D data
- npointsint
number of points in the output data
- intrangefloat
integration range in along intdir in degree. data will be integrated from -intrange .. +intrange
- intdir{‘omega’, ‘2theta’}, optional
integration direction: ‘omega’: sample rocking angle (default), ‘2theta’: scattering angle
- wlfloat or str, optional
wavelength used to determine angular integration positions
- Note:
Although applicable for any set of data, the extraction only makes sense when the data are aligned into the y/z-plane.
- Returns:
- tt, omttintndarray
omega-2theta scan coordinates (2theta values) and intensities
- used_maskndarray
mask of used data, shape is the same as the input intensity: True for points which contributed, False for all others
Examples
>>> ttcut, omtt_int, mask = get_radial_scan([qy, qz], inten, [2.0, 5.0], ... 250, intrange=0.1)
- xrayutilities.analysis.line_cuts.get_ttheta_scan(qpos, intensity, cutpos, npoints, intrange, **kwargs)[source]¶
extracts a 2theta scan from reciprocal space map data with integration along either the omega or radial direction. The coplanar diffraction geometry with qy and qz as diffraction plane is assumed. This is consistent with the coplanar geometry implemented in the HXRD-experiment class.
This function works for 2D and 3D input data in the same way!
- Parameters:
- qposlist of array-like objects
arrays of y, z (list with two components) or x, y, z (list with three components) momentum transfers
- intensityarray-like
2D or 3D array of reciprocal space intensity with shape equal to the qpos entries
- cutpostuple or list
y/z-position or x/y/z-position at which the line scan should be extracted. this must be have two entries for 2D data (z-position) and a three for 3D data
- npointsint
number of points in the output data
- intrangefloat
integration range in along intdir in degree. data will be integrated from -intrange .. +intrange
- intdir{‘omega’, ‘radial’}, optional
integration direction: ‘omega’: sample rocking angle (default), ‘radial’: omega-2theta direction
- wlfloat or str, optional
wavelength used to determine angular integration positions
- Note:
Although applicable for any set of data, the extraction only makes sense when the data are aligned into the y/z-plane.
- Returns:
- tt, ttintndarray
2theta scan coordinates and intensities
- used_maskndarray
mask of used data, shape is the same as the input intensity: True for points which contributed, False for all others
Examples
>>> ttcut, tt_int, mask = get_ttheta_scan([qy, qz], inten, [2.0, 5.0], ... 250, intrange=0.1)
xrayutilities.analysis.misc module¶
miscellaneous functions helpful in the analysis and experiment
- xrayutilities.analysis.misc.coplanar_intensity(mat, exp, hkl, thickness, thMono, sample_width=10, beam_width=1)[source]¶
Calculates the expected intensity of a Bragg peak from an epitaxial thin film measured in coplanar geometry (integration over omega and 2theta in angular space!)
- Parameters:
- matCrystal
Crystal instance for structure factor calculation
- expExperiment
Experimental(HXRD) class for the angle calculation
- hkllist, tuple or array-like
Miller indices of the peak to calculate
- thicknessfloat
film thickness in nm
- thMonofloat
Bragg angle of the monochromator (deg)
- sample_widthfloat, optional
width of the sample along the beam
- beam_widthfloat, optional
width of the beam in the same units as the sample size
- Returns:
- float
intensity of the peak
- xrayutilities.analysis.misc.getangles(peak, sur, inp)[source]¶
calculates the chi and phi angles for a given peak
- Parameters:
- peaklist or array-like
hkl for the peak of interest
- surlist or array-like
hkl of the surface
- inplist or array-like
inplane reference peak or direction
- Returns:
- list
[chi, phi] for the given peak on surface sur with inplane direction inp as reference
Examples
To get the angles for the -224 peak on a 111 surface type
>>> [chi, phi] = getangles([-2, 2, 4], [1, 1, 1], [2, 2, 4])
- xrayutilities.analysis.misc.getunitvector(chi, phi, ndir=(0, 0, 1), idir=(1, 0, 0))[source]¶
return unit vector determined by spherical angles and definition of the polar axis and inplane reference direction (phi=0)
- Parameters:
- chi, phifloat
spherical angles (polar and azimuthal) in degree
- ndirtuple, list or array-like
polar/z-axis (determines chi=0)
- idirtuple, list or array-like
azimuthal axis (determines phi=0)
xrayutilities.analysis.sample_align module¶
functions to help with experimental alignment during experiments, especially for experiments with linear and area detectors
- xrayutilities.analysis.sample_align.area_detector_calib(angle1, angle2, ccdimages, detaxis, r_i, plot=True, cut_off=0.7, start=(None, None, 1, 0, 0, 0, 0), fix=(False, False, True, False, False, False, False), fig=None, wl=None, plotlog=False, nwindow=50, debug=False)[source]¶
function to calibrate the detector parameters of an area detector it determines the detector tilt possible rotations and offsets in the detector arm angles
- Parameters:
- angle1array-like
outer detector arm angle
- angle2array-like
inner detector arm angle
- ccdimagesarray-like
images of the ccd taken at the angles given above
- detaxislist of str
detector arm rotation axis; default: [‘z+’, ‘y-‘]
- r_istr
primary beam direction [xyz][+-]; default ‘x+’
- plotbool, optional
flag to determine if results and intermediate results should be plotted; default: True
- cut_offfloat, optional
cut off intensity to decide if image is used for the determination or not; default: 0.7 = 70%
- starttuple, optional
sequence of start values of the fit for parameters, which can not be estimated automatically or might want to be fixed. These are: pwidth1, pwidth2, distance, tiltazimuth, tilt, detector_rotation, outerangle_offset. By default (None, None, 1, 0, 0, 0, 0) is used.
- fixtuple of bool
fix parameters of start (default: (False, False, True, False, False, False, False)) It is strongly recommended to either fix the distance or the pwidth1, 2 values.
- figFigure, optional
matplotlib figure used for plotting the error default: None (creates own figure)
- wlfloat or str
wavelength of the experiment in angstrom (default: config.WAVELENGTH) value does not really matter here but does affect the scaling of the error
- plotlogbool
flag to specify if the created error plot should be on log-scale
- nwindowint
window size for determination of the center of mass position after the center of mass of every full image is determined, the center of mass is determined again using a window of size nwindow in order to reduce the effect of hot pixels.
- debugbool
flag to specify that you want to see verbose output and saving of images to show if the CEN determination works
- xrayutilities.analysis.sample_align.area_detector_calib_hkl(sampleang, angle1, angle2, ccdimages, hkls, experiment, material, detaxis, r_i, plot=True, cut_off=0.7, start=(None, None, 1, 0, 0, 0, 0, 0, 0, 'config'), fix=(False, False, True, False, False, False, False, False, False, False), fig=None, plotlog=False, nwindow=50, debug=False)[source]¶
function to calibrate the detector parameters of an area detector it determines the detector tilt possible rotations and offsets in the detector arm angles
in this variant not only scans through the primary beam but also scans at a set of symmetric reflections can be used for the detector parameter determination. for this not only the detector parameters but in addition the sample orientation and wavelength need to be fit. Both images from the primary beam hkl = (0, 0, 0) and from a symmetric reflection hkl = (h, k, l) need to be given for a successful run.
- Parameters:
- sampleangarray-like
sample rocking angle (needed to align the reflections (same rotation direction as inner detector rotation)) other sample angle are not allowed to be changed during the scans
- angle1array-like
outer detector arm angle
- angle2array-like
inner detector arm angle
- ccdimagesarray-like
images of the ccd taken at the angles given above
- hklslist or array-like
hkl values for every image
- experimentExperiment
Experiment class object needed to get the UB matrix for the hkl peak treatment
- materialCrystal
material used as reference crystal
- detaxislist of str
detector arm rotation axis; default: [‘z+’, ‘y-‘]
- r_istr
primary beam direction [xyz][+-]; default ‘x+’
- plotbool, optional
flag to determine if results and intermediate results should be plotted; default: True
- cut_offfloat, optional
cut off intensity to decide if image is used for the determination or not; default: 0.7 = 70%
- starttuple, optional
sequence of start values of the fit for parameters, which can not be estimated automatically or might want to be fixed. These are: pwidth1, pwidth2, distance, tiltazimuth, tilt, detector_rotation, outerangle_offset, sampletilt, sampletiltazimuth, wavelength. By default (None, None, 1, 0, 0, 0, 0, 0, 0, ‘config’).
- fixtuple of bool
fix parameters of start (default: (False, False, True, False, False, False, False, False, False, False)) It is strongly recommended to either fix the distance or the pwidth1, 2 values.
- figFigure, optional
matplotlib figure used for plotting the error default: None (creates own figure)
- plotlogbool
flag to specify if the created error plot should be on log-scale
- nwindowint
window size for determination of the center of mass position after the center of mass of every full image is determined, the center of mass is determined again using a window of size nwindow in order to reduce the effect of hot pixels.
- debugbool
flag to specify that you want to see verbose output and saving of images to show if the CEN determination works
- xrayutilities.analysis.sample_align.fit_bragg_peak(om, tt, psd, omalign, ttalign, exphxrd, frange=(0.03, 0.03), peaktype='Gauss', plot=True)[source]¶
helper function to determine the Bragg peak position in a reciprocal space map used to obtain the position needed for correction of the data. the determination is done by fitting a two dimensional Gaussian (xrayutilities.math.Gauss2d) or Lorentzian (xrayutilities.math.Lorentz2d)
PLEASE ALWAYS CHECK THE RESULT CAREFULLY!
- Parameters:
- om, ttarray-like
angular coordinates of the measurement either with size of psd or of psd.shape[0]
- psdarray-like
intensity values needed for fitting
- omalignfloat
aligned omega value, used as first guess in the fit
- ttalignfloat
aligned two theta values used as first guess in the fit these values are also used to set the range for the fit: the peak should be within +/-frangeAA^{-1} of those values
- exphxrdExperiment
experiment class used for the conversion between angular and reciprocal space.
- frangetuple of float, optional
data range used for the fit in both directions (see above for details default:(0.03, 0.03) unit: AA^{-1})
- peaktype{‘Gauss’, ‘Lorentz’}
peak type to fit
- plotbool, optional
if True (default) function will plot the result of the fit in comparison with the measurement.
- Returns:
- omfit, ttfitfloat
fitted angular values
- paramslist
fit parameters (of the Gaussian/Lorentzian)
- covariancendarray
covariance matrix of the fit parameters
- xrayutilities.analysis.sample_align.linear_detector_calib(angle, mca_spectra, **keyargs)[source]¶
function to calibrate the detector distance/channel per degrees for a straight linear detector mounted on a detector arm
- Parameters:
- anglearray-like
array of angles in degree of measured detector spectra
- mca_spectraarray-like
corresponding detector spectra (shape: (len(angle), Nchannels)
- r_istr, optional
primary beam direction as vector [xyz][+-]; default: ‘y+’
- detaxisstr, optional
detector arm rotation axis [xyz][+-]; default: ‘x+’
- Returns:
- pixelwidthfloat
width of the pixel at one meter distance, pixelwidth is negative in case the hit channel number decreases upon an increase of the detector angle
- center_channelfloat
central channel of the detector
- detector_tiltfloat, optional
if usetilt=True the fitted tilt of the detector is also returned
Note
L/pixelwidth*pi/180 ~= channel/degree, with the sample detector distance L
- The function also prints out how a linear detector can be initialized using
- the results obtained from this calibration. Carefully check the results
- Other Parameters:
- plotbool
flag to specify if a visualization of the fit should be done
- usetiltbool
whether to use model considering a detector tilt, i.e. deviation angle of the pixel direction from orthogonal to the primary beam (default: True)
See also
psd_chdeg
low level function with more configurable options
- xrayutilities.analysis.sample_align.miscut_calc(phi, aomega, zeros=None, omega0=None, plot=True)[source]¶
function to calculate the miscut direction and miscut angle of a sample by fitting a sinusoidal function to the variation of the aligned omega values of more than two reflections. The function can also be used to fit reflectivity alignment values in various azimuths.
- Parameters:
- philist, tuple or array-like
azimuths in which the reflection was aligned (deg)
- aomegalist, tuple or array-like
aligned omega values (deg)
- zeroslist, tuple or array-like, optional
angles at which surface is parallel to the beam (deg). For the analysis the angles (aomega - zeros) are used.
- omega0float, optional
if specified the nominal value of the reflection is not included as fit parameter, but is fixed to the specified value. This value is MANDATORY if ONLY TWO AZIMUTHs are given.
- plotbool, optional
flag to specify if a visualization of the fit is wanted. default: True
- Returns:
- omega0float
the omega value of the reflection should be close to the nominal one
- phi0float
the azimuth in which the primary beam looks upstairs
- miscutfloat
amplitude of the sinusoidal variation == miscut angle
- xrayutilities.analysis.sample_align.psd_chdeg(angles, channels, stdev=None, usetilt=True, plot=True, datap='xk', modelline='--r', modeltilt='-b', fignum=None, mlabel='fit', mtiltlabel='fit w/tilt', dlabel='data', figtitle=True)[source]¶
function to determine the channels per degree using a linear fit of the function nchannel = center_ch+chdeg*tan(angles) or the equivalent including a detector tilt
- Parameters:
- anglesarray-like
detector angles for which the position of the beam was measured
- channelsarray-like
detector channels where the beam was found
- stdevarray-like, optional
standard deviation of the beam position
- plotbool, optional
flag to specify if a visualization of the fit should be done
- usetiltbool, optional
whether to use model considering a detector tilt, i.e. deviation angle of the pixel direction from orthogonal to the primary beam (default: True)
- Returns:
- pixelwidthfloat
the width of one detector channel @ 1m distance, which is negative in case the hit channel number decreases upon an increase of the detector angle.
- centerchfloat
center channel of the detector
- tiltfloat
tilt of the detector from perpendicular to the beam (will be zero in case of usetilt=False)
Note
L/pixelwidth*pi/180 = channel/degree for large detector distance with the sample detector disctance L
- Other Parameters:
- datapstr, optional
plot format of data points
- modellinestr, optional
plot format of modelline
- modeltiltstr, optional
plot format of modeltilt
- fignumint or str, optional
figure number to use for the plot
- mlabelstr
label of the model w/o tilt to be used in the plot
- mtiltlabelstr
label of the model with tilt to be used in the plot
- dlabelstr
label of the data line to be used in the plot
- figtitlebool
flag to tell if the figure title should show the fit parameters
- xrayutilities.analysis.sample_align.psd_refl_align(primarybeam, angles, channels, plot=True)[source]¶
function which calculates the angle at which the sample is parallel to the beam from various angles and detector channels from the reflected beam. The function can be used during the half beam alignment with a linear detector.
- Parameters:
- primarybeamint
primary beam channel number
- angleslist or array-like
incidence angles
- channelslist or array-like
corresponding detector channels
- plotbool, optional
flag to specify if a visualization of the fit is wanted default : True
- Returns:
- float
angle at which the sample is parallel to the beam
Examples
>>> zeroangle = psd_refl_align(500, [0, 0.1, 0.2, 0.3], ... [550, 600, 640, 700]) XU.analysis.psd_refl_align: sample is parallel to beam at goniometer angle -0.0986 (R^2=0.9942)
Module contents¶
xrayutilities.analysis is a package for assisting with the analysis of x-ray diffraction data, mainly reciprocal space maps
Routines for obtaining line cuts from gridded reciprocal space maps are offered, with the ability to integrate the intensity perpendicular to the line cut direction.