pyFAI.geometry package
Module contents
This sub-package contains mostly a large Geometry class (defined in core) which is in charge of:
calculating the geometry, i.e. the position in the detector space of each pixel of the detector
manages caches to store intermediate results
There are several other helper classes which provide conversion to other software
NOTA: The Geometry class is not a “transformation class” which would take a detector and transform it. It is rather a description of the experimental setup.
pyFAI.geometry.core module
This modules contains only one (large) class in charge of:
calculating the geometry, i.e. the position in the detector space of each pixel of the detector
manages caches to store intermediate results
NOTA: The Geometry class is not a “transformation class” which would take a detector and transform it. It is rather a description of the experimental setup.
- class pyFAI.geometry.core.Geometry(dist=1, poni1=0, poni2=0, rot1=0, rot2=0, rot3=0, pixel1=None, pixel2=None, splineFile=None, detector=None, wavelength=None, orientation=0)
Bases:
object
This class is the parent-class of azimuthal integrator.
This class contains a detector (using composition) which provides the position of all pixels, or only a limited set of pixel indices. The Geometry class is responsible for translating/rotating those pixel to their position in reference to the sample position. The description of the experimental setup is inspired by the work of P. Boesecke
Detector is assumed to be corrected from “raster orientation” effect. It is not addressed here but rather in the Detector object or at read time. Considering there is no tilt:
Detector fast dimension (dim2) is supposed to be horizontal (dimension X of the image)
Detector slow dimension (dim1) is supposed to be vertical, upwards (dimension Y of the image)
The third dimension is chose such as the referential is orthonormal, so dim3 is along incoming X-ray beam
Demonstration of the equation done using Mathematica:
- __init__(dist=1, poni1=0, poni2=0, rot1=0, rot2=0, rot3=0, pixel1=None, pixel2=None, splineFile=None, detector=None, wavelength=None, orientation=0)
- Parameters:
dist – distance sample - detector plan (orthogonal distance, not along the beam), in meter.
poni1 – coordinate of the point of normal incidence along the detector’s first dimension, in meter
poni2 – coordinate of the point of normal incidence along the detector’s second dimension, in meter
rot1 – first rotation from sample ref to detector’s ref, in radians
rot2 – second rotation from sample ref to detector’s ref, in radians
rot3 – third rotation from sample ref to detector’s ref, in radians
pixel1 (float) – Deprecated. Pixel size of the fist dimension of the detector, in meter. If both pixel1 and pixel2 are not None, detector pixel size is overwritten. Prefer defining the detector pixel size on the provided detector object. Prefer defining the detector pixel size on the provided detector object (
detector.pixel1 = 5e-6
).pixel2 (float) – Deprecated. Pixel size of the second dimension of the detector, in meter. If both pixel1 and pixel2 are not None, detector pixel size is overwritten. Prefer defining the detector pixel size on the provided detector object (
detector.pixel2 = 5e-6
).splineFile (str) – Deprecated. File containing the geometric distortion of the detector. If not None, pixel1 and pixel2 are ignored and detector spline is overwritten. Prefer defining the detector spline manually (
detector.splineFile = "file.spline"
).detector (str or pyFAI.Detector) – name of the detector or Detector instance. String description is deprecated. Prefer using the result of the detector factory:
pyFAI.detector_factory("eiger4m")
wavelength (float) – Wave length used in meter
orientation (int) – orientation of the detector, see pyFAI.detectors.orientation.Orientation
- array_from_unit(shape=None, typ='center', unit=2th_deg, scale=True)
Generate an array of position in different dimentions (R, Q, 2Theta …)
- Parameters:
shape (ndarray.shape) – shape of the expected array, leave it to None for safety
typ (str) – “center”, “corner” or “delta”
unit (pyFAI.units.Unit or 2-tuple of them (valid only for corner coordinates calculation) – can be any valid unit (found in units.AZIMUTHAL_UNITS or units.RADIAL_UNITS)
scale – set to False for returning the internal representation (S.I. often) which is faster
- Returns:
R, Q or 2Theta array depending on unit
- Return type:
ndarray
- calc_pos_zyx(d0=None, d1=None, d2=None, param=None, corners=False, use_cython=True, do_parallax=False)
Calculate the position of a set of points in space in the sample’s centers referential.
This is usually used for calculating the pixel position in space.
Nota: dim3 is the same as dim0
- Parameters:
d0 – altitude on the point compared to the detector (i.e. z), may be None
d1 – position on the detector along the slow dimension (i.e. y)
d2 – position on the detector along the fastest dimension (i.e. x)
corners – return positions on the corners (instead of center)
use_cython – set to False to validate using pure numpy
do_parallax – position should be corrected for parallax effect
- Returns:
3-tuple of nd-array, with dim0=along the beam, dim1=along slowest dimension dim2=along fastest dimension
- calc_transmission(t0, shape=None)
Defines the absorption correction for a phosphor screen or a scintillator from t0, the normal transmission of the screen.
\[ \begin{align}\begin{aligned}Icor = \frac{Iobs(1-t0)}{1-exp(ln(t0)/cos(incidence))}\\let_t = \frac{1-exp(ln(t0)/cos(incidence))}{1 - t0}\end{aligned}\end{align} \]See reference on: J. Appl. Cryst. (2002). 35, 356–359 G. Wu et al. CCD phosphor
- Parameters:
t0 – value of the normal transmission (from 0 to 1)
shape – shape of the array
- Returns:
actual
- calcfrom1d(tth, I, shape=None, mask=None, dim1_unit=2th_deg, correctSolidAngle=True, dummy=0.0, polarization_factor=None, polarization_axis_offset=0, dark=None, flat=None)
Computes a 2D image from a 1D integrated profile
- Parameters:
tth – 1D array with radial unit, this array needs to be ordered
I – scattering intensity, corresponding intensity
shape – shape of the image (if not defined by the detector)
dim1_unit – unit for the “tth” array
correctSolidAngle –
dummy – value for masked pixels
polarization_factor – set to true to use previously used value
polarization_axis_offset – axis_offset to be send to the polarization method
dark – dark current correction
flat – flatfield corrction
- Returns:
2D image reconstructed
- calcfrom2d(I, tth, chi, shape=None, mask=None, dim1_unit=2th_deg, dim2_unit=chi_deg, correctSolidAngle=True, dummy=0.0, polarization_factor=None, polarization_axis_offset=0, dark=None, flat=None)
Computes a 2D image from a cake / 2D integrated image
- Parameters:
I – scattering intensity, as an image n_tth, n_chi
tth – 1D array with radial unit, this array needs to be ordered
chi – 1D array with azimuthal unit, this array needs to be ordered
shape – shape of the image (if not defined by the detector)
dim1_unit – unit for the “tth” array
dim2_unit – unit for the “chi” array
correctSolidAngle –
dummy – value for masked pixels
polarization_factor – set to true to use previously used value
polarization_axis_offset – axis_offset to be send to the polarization method
dark – dark current correction
flat – flatfield corrction
- Returns:
2D image reconstructed
- center_array(shape=None, unit='2th_deg', scale=True)
Generate a 2D array of the given shape with (i,j) (radial angle ) for all elements.
- Parameters:
shape (2-tuple of integer) – expected shape
unit – string like “2th_deg” or an instance of pyFAI.units.Unit
scale – set to False for returning the internal representation (S.I. often) which is faster.
- Returns:
2d array of given shape
- check_chi_disc(azimuth_range)
Check the position of the \(\chi\) discontinuity
- Parameters:
range – range of chi for the integration
- Returns:
True if there is a problem
- chi(d1, d2, path='cython')
Calculate the chi (azimuthal angle) for the centre of a pixel at coordinate d1, d2. Conversion to lab coordinate system is performed in calc_pos_zyx.
- Parameters:
d1 (float or array of them) – pixel coordinate along the 1st dimention (C convention)
d2 (float or array of them) – pixel coordinate along the 2nd dimention (C convention)
path – can be “tan” (i.e via numpy) or “cython”
- Returns:
chi, the azimuthal angle in rad
- chiArray(shape=None)
Generate an array of azimuthal angle chi(i,j) for all elements in the detector.
Azimuthal angles are in radians
Nota: Refers to the pixel centers !
- Parameters:
shape – the shape of the chi array
- Returns:
the chi array as numpy.ndarray
- chi_corner(d1, d2)
Calculate the chi (azimuthal angle) for the corner of a pixel at coordinate d1,d2 which in the lab ref has coordinate:
- Parameters:
d1 (float or array of them) – pixel coordinate along the 1st dimention (C convention)
d2 (float or array of them) – pixel coordinate along the 2nd dimention (C convention)
- Returns:
chi, the azimuthal angle in rad
- property chia
chi array in cache
- cornerArray(shape=None)
Generate a 4D array of the given shape with (i,j) (radial angle 2th, azimuthal angle chi ) for all elements.
- Parameters:
shape (2-tuple of integer) – expected shape
- Returns:
3d array with shape=(*shape,4,2) the two elements are: * dim3[0]: radial angle 2th * dim3[1]: azimuthal angle chi
- cornerQArray(shape=None)
Generate a 3D array of the given shape with (i,j) (azimuthal angle) for all elements.
- Parameters:
shape (2-tuple of integer) – expected shape
- Returns:
3d array with shape=(*shape,4,2) the two elements are (scattering vector q, azimuthal angle chi)
- cornerRArray(shape=None)
Generate a 3D array of the given shape with (i,j) (azimuthal angle) for all elements.
- Parameters:
shape (2-tuple of integer) – expected shape
- Returns:
3d array with shape=(*shape,4,2) the two elements are (radial distance, azimuthal angle chi)
- cornerRd2Array(shape=None)
Generate a 3D array of the given shape with (i,j) (azimuthal angle) for all elements.
- Parameters:
shape (2-tuple of integer) – expected shape
- Returns:
3d array with shape=(*shape,4,2) the two elements are (reciprocal spacing squared, azimuthal angle chi)
- corner_array(shape=None, unit=None, use_cython=True, scale=True)
Generate a 3D array of the given shape with (i,j) (radial angle 2th, azimuthal angle chi ) for all elements.
- Parameters:
shape (2-tuple of integer) – expected shape
unit – string like “2th_deg” or an instance of pyFAI.units.Unit
use_cython – set to False to use the slower Python path (for tests)
scale – set to False for returning the internal representation (S.I. often) which is faster
- Returns:
3d array with shape=(*shape,4,2) the two elements are: - dim3[0]: radial angle 2th, q, r… - dim3[1]: azimuthal angle chi
- property correct_SA_spline
- cos_incidence(d1, d2, path='cython')
Calculate the cosinus of the incidence angle (alpha) for current pixels (P). The poni being the point of normal incidence, it’s incidence angle is \({\alpha} = 0\) hence \(cos({\alpha}) = 1\).
Works also for non-flat detectors where the normal of the pixel is considered.
- Parameters:
d1 – 1d or 2d set of points in pixel coord
d2 – 1d or 2d set of points in pixel coord
path – can be “cython”, “numexpr” or “numpy” (fallback).
- Returns:
cosine of the incidence angle
- del_chia()
- del_dssa()
- del_qa()
- del_ra()
- del_ttha()
- delta2Theta(shape=None)
Generate a 2D array of the given shape with the max distance between the center and any corner in 2 theta angle in radians
- Parameters:
shape – The shape of the detector array: 2-tuple of integer
- Returns:
2D-array containing the max delta angle between a pixel center and any corner in 2theta-angle (rad)
- deltaChi(shape=None, use_cython=True)
Generate a 2D array of the given shape with the max distance between the center and any corner in chi-angle (rad)
- Parameters:
shape – The shape of the detector array: 2-tuple of integer
- Returns:
2D-array containing the max delta angle between a pixel center and any corner in chi-angle (rad)
- deltaQ(shape=None)
Generate a 2D array of the given shape with the max distance between the center and any corner in q_vector unit (nm^-1)
- Parameters:
shape – The shape of the detector array: 2-tuple of integer
- Returns:
array 2D containing the max delta Q between a pixel center and any corner in q_vector unit (nm^-1)
- deltaR(shape=None)
Generate a 2D array of the given shape with (i,j) with the max distance between the center and any corner in radius unit (m)
- Parameters:
shape – The shape of the detector array: 2-tuple of integer
- Returns:
array 2D containing the max delta Q between a pixel center and any corner in distance (m)
- deltaRd2(shape=None)
Generate a 2D array of the given shape with the max distance between the center and any corner in unit: reciprocal spacing squarred (1/nm^2)
- Parameters:
shape – The shape of the detector array: 2-tuple of integer
- Returns:
array 2D containing the max delta (d*)^2 between a pixel center and any corner in reciprocal spacing squarred (1/nm^2)
- delta_array(shape=None, unit='2th_deg', scale=False)
Generate a 2D array of the given shape with the delta-radius for all elements.
- Parameters:
shape (2-tuple of integer) – expected shape
unit – string like “2th_deg” or an instance of pyFAI.units.Unit
scale – set to False for returning the internal representation (S.I. often) which is faster
- Returns:
2D array
- diffSolidAngle(d1, d2)
Calculate the solid angle of the current pixels (P) versus the PONI (C)
\[ \begin{align}\begin{aligned}dOmega = \frac{Omega(P)}{Omega(C)} = \frac{A \cdot cos(a)}{SP^2} \cdot \frac{SC^2}{A \cdot cos(0)} = \frac{3}{cos(a)} = \frac{SC^3}{SP^3}\\cos(a) = \frac{SC}{SP}\end{aligned}\end{align} \]- Parameters:
d1 – 1d or 2d set of points
d2 – 1d or 2d set of points (same size&shape as d1)
- Returns:
solid angle correction array
- property dist
- property dssa
solid angle array in cache
- property energy
- getCXI()
Export the geometry in the CXI format as defined in https://raw.githubusercontent.com/cxidb/CXI/master/cxi_file_format.pdf
- Returns:
dict with the structure of a CXI file to be written into HDF5
- getFit2D()
Export geometry setup with the geometry of Fit2D
- Returns:
dict with parameters compatible with fit2D geometry
- getImageD11()
Export the current geometry in ImageD11 format. Please refer to the documentation in doc/source/geometry_conversion.rst for the orientation and units of those values.
- Returns:
an Ordered dict with those parameters: distance 294662.658 #in nm o11 1 o12 0 o21 0 o22 -1 tilt_x 0.00000 tilt_y -0.013173 tilt_z 0.002378 wavelength 0.154 y_center 1016.328171 y_size 48.0815 z_center 984.924425 z_size 46.77648
- getPyFAI()
Export geometry setup with the geometry of PyFAI
Deprecated, please use get/set_config which is cleaner! when it comes to detector specification
- Returns:
dict with the parameter-set of the PyFAI geometry
- getSPD()
get the SPD like parameter set: For geometry description see Peter Boesecke J.Appl.Cryst.(2007).40, s423–s427
Basically the main difference with pyFAI is the order of the axis which are flipped
- Returns:
dictionnary with those parameters: SampleDistance: distance from sample to detector at the PONI (orthogonal projection) Center_1, pixel position of the PONI along fastest axis Center_2: pixel position of the PONI along slowest axis Rot_1: rotation around the fastest axis (x) Rot_2: rotation around the slowest axis (y) Rot_3: rotation around the axis ORTHOGONAL to the detector plan PSize_1: pixel size in meter along the fastest dimention PSize_2: pixel size in meter along the slowst dimention splineFile: name of the file containing the spline BSize_1: pixel binning factor along the fastest dimention BSize_2: pixel binning factor along the slowst dimention WaveLength: wavelength used in meter
- get_chia()
- get_config()
return the configuration as a dictionnary
- Returns:
dictionary with the current configuration
- get_correct_solid_angle_for_spline()
- get_dist()
- get_dssa()
- get_energy()
- get_mask()
- get_maskfile()
- get_parallax()
- get_pixel1()
- get_pixel2()
- get_poni1()
- get_poni2()
- get_qa()
- get_ra()
- get_rot1()
- get_rot2()
- get_rot3()
- get_shape(shape=None)
Guess what is the best shape ….
- Parameters:
shape – force this value (2-tuple of int)
- Returns:
2-tuple of int
- get_spline()
- get_splineFile()
- get_ttha()
- get_wavelength()
- guess_npt_rad()
calculate the number of pixels from the beam-center to the corner the further away from it.
- Returns:
this distance as a number of pixels.
It is a good guess of the number of bins to be used without oversampling too much the data for azimuthal integration
- load(filename)
Load the refined parameters from a file.
- Parameters:
filename (string) – name of the file to load
- Returns:
itself with updated parameters
- make_headers(type_='list')
Create a configuration / header for the integrated data
- Parameters:
type – can be “list” or “dict”
- Returns:
the header with the proper format
- property mask
- property maskfile
- normalize_azimuth_range(azimuth_range)
Convert the azimuth range from degrees to radians
This method takes care of the position of the discontinuity and adapts the range accordingly!
- Parameters:
azimuth_range – 2-tuple of float in degrees
- Returns:
2-tuple of float in radians in a range such to avoid the discontinuity
- oversampleArray(myarray)
- property parallax
- property pixel1
- property pixel2
- polarization(shape=None, factor=None, axis_offset=0, with_checksum=False, path='numexpr')
Calculate the polarization correction accoding to the polarization factor:
- If the polarization factor is None,
the correction is not applied (returns 1)
- If the polarization factor is 0 (circular polarization),
the correction correspond to (1+(cos2θ)^2)/2
- If the polarization factor is 1 (linear horizontal polarization),
there is no correction in the vertical plane and a node at 2th=90, chi=0
- If the polarization factor is -1 (linear vertical polarization),
there is no correction in the horizontal plane and a node at 2th=90, chi=90
If the polarization is elliptical, the polarization factor varies between -1 and +1.
The axis_offset parameter allows correction for the misalignement of the polarization plane (or ellipse main axis) and the the detector’s X axis.
- Parameters:
shape – the shape of the array, can be guessed most of the time from the detector definition
factor – (Ih-Iv)/(Ih+Iv): varies between 0 (circular/random polarization) and 1 (where division by 0 could occure at 2th=90, chi=0)
axis_offset – Angle between the polarization main axis and detector’s X direction (in radians !!!)
with_checksum – calculate also the checksum (used with OpenCL integration)
path – set to numpy to enforce the use of numpy, else uses numexpr (mutithreaded)
- Returns:
2D array with polarization correction (normalization) array or namedtuple if with_checksum
- property poni1
- property poni2
- positionArray(*arg, **kwarg)
Deprecated version of
position_array()
, left for compatibility see doc of position_array
- position_array(shape=None, corners=False, dtype=<class 'numpy.float64'>, use_cython=True, do_parallax=False)
Generate an array for the pixel position given the shape of the detector.
if corners is False, the coordinates of the center of the pixel is returned in an array of shape: (shape[0], shape[1], 3) where the 3 coordinates are: * z: along incident beam, * y: to the top/sky, * x: towards the center of the ring
If is True, the corner of each pixels are then returned. the output shape is then (shape[0], shape[1], 4, 3)
- Parameters:
shape – shape of the array expected
corners – set to true to receive a (…,4,3) array of corner positions
dtype – output format requested. Double precision is needed for fitting the geometry
use_cython ((bool)) – set to false to test the Python path (slower)
do_parallax – correct for parallax effect (if parametrized)
- Returns:
3D coodinates as nd-array of size (…,3) or (…,3) (default)
Nota: this value is not cached and actually generated on demand (costly)
- qArray(shape=None)
Generate an array of the given shape with q(i,j) for all elements.
- qCornerFunct(d1, d2)
Calculate the q_vector for any pixel corner (in nm^-1)
- Parameters:
shape – expected shape of the detector
- qFunction(d1, d2, param=None, path='cython')
Calculates the q value for the center of a given pixel (or set of pixels) in nm-1
q = 4pi/lambda sin( 2theta / 2 )
- Parameters:
d1 (scalar or array of scalar) – position(s) in pixel in first dimension (c order)
d2 (scalar or array of scalar) – position(s) in pixel in second dimension (c order)
- Returns:
q in in nm^(-1)
- Return type:
float or array of floats.
- property qa
Q array in cache
- quaternion(param=None)
Calculate the quaternion associated to the current rotations from rot1, rot2, rot3.
Uses the transformations-library from C. Gohlke
- Parameters:
param – use this set of parameters instead of the default one.
- Returns:
numpy array with 4 elements [w, x, y, z]
- rArray(shape=None)
Generate an array of the given shape with r(i,j) for all elements; The radius r being in meters.
- Parameters:
shape – expected shape of the detector
- Returns:
2d array of the given shape with radius in m from beam center on detector.
- rCornerFunct(d1, d2)
Calculate the radius array for any pixel corner (in m)
- rFunction(d1, d2, param=None, path='cython')
Calculates the radius value for the center of a given pixel (or set of pixels) in m
r = distance to the incident beam
- Parameters:
d1 (scalar or array of scalar) – position(s) in pixel in first dimension (c order)
d2 (scalar or array of scalar) – position(s) in pixel in second dimension (c order)
- Returns:
r in in m
- Return type:
float or array of floats.
- property ra
R array in cache
- rd2Array(shape=None)
Generate an array of the given shape with (d*(i,j))^2 for all pixels.
d*^2 is the reciprocal spacing squared in inverse nm squared
- Parameters:
shape – expected shape of the detector
- Returns:
2d array of the given shape with reciprocal spacing squared
- read(filename)
Load the refined parameters from a file.
- Parameters:
filename (string) – name of the file to load
- Returns:
itself with updated parameters
- reset(collect_garbage=True)
reset most arrays that are cached: used when a parameter changes.
- Parameters:
collect_garbage – set to False to prevent garbage collection, faster
- property rot1
- property rot2
- property rot3
- rotation_matrix(param=None)
Compute and return the detector tilts as a single rotation matrix
Corresponds to rotations about axes 1 then 2 then 3 (=> 0 later on) For those using spd (PB = Peter Boesecke), tilts relate to this system (JK = Jerome Kieffer) as follows: JK1 = PB2 (Y) JK2 = PB1 (X) JK3 = PB3 (Z) …slight differences will result from the order FIXME: make backwards and forwards converter helper function
axis1 is vertical and perpendicular to beam axis2 is horizontal and perpendicular to beam axis3 is along the beam, becomes axis0 see: http://pyfai.readthedocs.io/en/latest/geometry.html#detector-position or ../doc/source/img/PONI.png
- Parameters:
param (list of float) – list of geometry parameters, defaults to self.param uses elements [3],[4],[5]
- Returns:
rotation matrix
- Return type:
3x3 float array
- save(filename)
Save the geometry parameters.
- Parameters:
filename (string) – name of the file where to save the parameters
- setCXI(dico)
Set the geometry of the azimuthal integrator from a CXI data structure (as dict)
- Parameters:
dico – dictionary with CXI information
- setChiDiscAtPi()
Set the position of the discontinuity of the chi axis between -pi and +pi. This is the default behavour
- setChiDiscAtZero()
Set the position of the discontinuity of the chi axis between 0 and 2pi. By default it is between pi and -pi
- setFit2D(directDist, centerX, centerY, tilt=0.0, tiltPlanRotation=0.0, pixelX=None, pixelY=None, splineFile=None, detector=None, wavelength=None)
Set the Fit2D-like parameter set: For geometry description see HPR 1996 (14) pp-240 https://doi.org/10.1080/08957959608201408
Warning: Fit2D flips automatically images depending on their file-format. By reverse engineering we noticed this behavour for Tiff and Mar345 images (at least). To obtaine correct result you will have to flip images using numpy.flipud.
- Parameters:
direct – direct distance from sample to detector along the incident beam (in millimeter as in fit2d)
tilt – tilt in degrees
tiltPlanRotation – Rotation (in degrees) of the tilt plan arround the Z-detector axis * 0deg -> Y does not move, +X goes to Z<0 * 90deg -> X does not move, +Y goes to Z<0 * 180deg -> Y does not move, +X goes to Z>0 * 270deg -> X does not move, +Y goes to Z>0
pixelX,pixelY – as in fit2d they ar given in micron, not in meter
centerY (centerX,) – pixel position of the beam center
splineFile – name of the file containing the spline
detector – name of the detector or detector object
- setImageD11(param)
Set the geometry from the parameter set which contains distance, o11, o12, o21, o22, tilt_x, tilt_y tilt_z, wavelength, y_center, y_size, z_center and z_size. Please refer to the documentation in doc/source/geometry_conversion.rst for the orientation and units of those values.
- Parameters:
param – dict with the values to set.
- setOversampling(iOversampling)
set the oversampling factor
- setPyFAI(**kwargs)
set the geometry from a pyFAI-like dict
Deprecated, please use get/set_config which is cleaner! when it comes to detector specification
- setSPD(SampleDistance, Center_1, Center_2, Rot_1=0, Rot_2=0, Rot_3=0, PSize_1=None, PSize_2=None, splineFile=None, BSize_1=1, BSize_2=1, WaveLength=None)
Set the SPD like parameter set: For geometry description see Peter Boesecke J.Appl.Cryst.(2007).40, s423–s427
Basically the main difference with pyFAI is the order of the axis which are flipped
- Parameters:
SampleDistance – distance from sample to detector at the PONI (orthogonal projection)
Center_1 – pixel position of the PONI along fastest axis
Center_2 – pixel position of the PONI along slowest axis
Rot_1 – rotation around the fastest axis (x)
Rot_2 – rotation around the slowest axis (y)
Rot_3 – rotation around the axis ORTHOGONAL to the detector plan
PSize_1 – pixel size in meter along the fastest dimention
PSize_2 – pixel size in meter along the slowst dimention
splineFile – name of the file containing the spline
BSize_1 – pixel binning factor along the fastest dimention
BSize_2 – pixel binning factor along the slowst dimention
WaveLength – wavelength used
- set_chia(_)
- set_config(config)
Set the config of the geometry and of the underlying detector
- Parameters:
config – dictionary with the configuration
- Returns:
itself
- set_correct_solid_angle_for_spline(value)
- set_dist(value)
- set_dssa(_)
- set_energy(energy)
Set the energy in keV
- set_mask(mask)
- set_maskfile(maskfile)
- set_parallax(value)
- set_param(param)
set the geometry from a 6-tuple with dist, poni1, poni2, rot1, rot2, rot3
- set_pixel1(pixel1)
- set_pixel2(pixel2)
- set_poni1(value)
- set_poni2(value)
- set_qa(_)
- set_ra(_)
- set_rot1(value)
- set_rot2(value)
- set_rot3(value)
- set_rot_from_quaternion(w, x, y, z)
Quaternions are convieniant ways to represent 3D rotation This method allows to define rot1(left-handed), rot2(left-handed) and rot3 (right handed) as definied in the documentation from a quaternion, expressed in the right handed (x1, x2, x3) basis set.
Uses the transformations-library from C. Gohlke
- Parameters:
w – Real part of the quaternion (correspond to cos alpha/2)
x – Imaginary part of the quaternion, correspond to u1*sin(alpha/2)
y – Imaginary part of the quaternion, correspond to u2*sin(alpha/2)
z – Imaginary part of the quaternion, correspond to u3*sin(alpha/2)
- set_spline(spline)
- set_splineFile(splineFile)
- set_ttha(_)
- set_wavelength(value)
Set the wavelength in meter!
- sin_incidence(d1, d2, path='cython')
Calculate the sinus of the incidence angle (alpha) for current pixels (P). The poni being the point of normal incidence, it’s incidence angle is \({\alpha} = 0\) hence \(sin({\alpha}) = 0\).
Works also for non-flat detectors where the normal of the pixel is considered.
- Parameters:
d1 – 1d or 2d set of points in pixel coord
d2 – 1d or 2d set of points in pixel coord
path – can be “cython”, “numexpr” or “numpy” (fallback).
- Returns:
cosine of the incidence angle
- classmethod sload(filename)
A static method combining the constructor and the loader from a file
- Parameters:
filename (string) – name of the file to load
- Returns:
instance of Geometry of AzimuthalIntegrator set-up with the parameter from the file.
- solidAngleArray(shape=None, order=3, absolute=False, with_checksum=False)
Generate an array for the solid angle correction given the shape of the detector.
solid_angle = cos(incidence)^3
- Parameters:
shape – shape of the array expected
order – should be 3, power of the formula just obove
absolute – the absolute solid angle is calculated as:
with_checksum – returns the array _and_ its checksum as a 2-tuple
SA = pix1*pix2/dist^2 * cos(incidence)^3
- property spline
- property splineFile
- tth(d1, d2, param=None, path='cython')
Calculates the 2theta value for the center of a given pixel (or set of pixels)
- Parameters:
d1 (scalar or array of scalar) – position(s) in pixel in first dimension (c order)
d2 (scalar or array of scalar) – position(s) in pixel in second dimension (c order)
path – can be “cos”, “tan” or “cython”
- Returns:
2theta in radians
- Return type:
float or array of floats.
- tth_corner(d1, d2)
Calculates the 2theta value for the corner of a given pixel (or set of pixels)
- Parameters:
d1 (scalar or array of scalar) – position(s) in pixel in first dimension (c order)
d2 (scalar or array of scalar) – position(s) in pixel in second dimension (c order)
- Returns:
2theta in radians
- Return type:
floar or array of floats.
- property ttha
2theta array in cache
- twoThetaArray(shape=None)
Generate an array of two-theta(i,j) in radians for each pixel in detector
the 2theta array values are in radians
- Parameters:
shape – shape of the detector
- Returns:
array of 2theta position in radians
- property wavelength
- write(filename)
Save the geometry parameters.
- Parameters:
filename (string) – name of the file where to save the parameters
pyFAI.geometry.cxi module
This modules contains helper function to convert to/from CXI format described in:
https://raw.githubusercontent.com/cxidb/CXI/master/cxi_file_format.pdf
ToDo
pyFAI.geometry.fit2d module
This modules contains helper function to convert to/from FIT2D geometry
- class pyFAI.geometry.fit2d.Fit2dGeometry(directDist: float = None, centerX: float = None, centerY: float = None, tilt: float = 0.0, tiltPlanRotation: float = 0.0, pixelX: float = None, pixelY: float = None, splineFile: str = None, detector: Detector = None, wavelength: float = None)
Bases:
NamedTuple
This object represents the geometry as configured in Fit2D
- Parameters:
directDist – Distance from sample to the detector along the incident beam in mm. The detector may be extrapolated when tilted.
centerX – Position of the beam-center on the detector in pixels, along the fastest axis of the image.
centerY – Position of the beam-center on the detector in pixels, along the slowest axis of the image.
tilt – Angle of tilt of the detector in degrees
tiltPlanRotation – Direction of the tilt (unefined when tilt is 0)
detector – Detector definition as is pyFAI.
wavelength – Wavelength of the beam in Angstrom
- centerX: float
Alias for field number 1
- centerY: float
Alias for field number 2
- directDist: float
Alias for field number 0
- pixelX: float
Alias for field number 5
- pixelY: float
Alias for field number 6
- splineFile: str
Alias for field number 7
- tilt: float
Alias for field number 3
- tiltPlanRotation: float
Alias for field number 4
- wavelength: float
Alias for field number 9
- pyFAI.geometry.fit2d.convert_from_Fit2d(f2d)
Import the geometry from Fit2D
- Parameters:
f2d – instance of Fit2dGeometry
- Returns:
PoniFile instance
- pyFAI.geometry.fit2d.convert_to_Fit2d(poni)
Convert a Geometry|PONI object to the geometry of Fit2D Please see the doc from Fit2dGeometry
- Parameters:
poni – azimuthal integrator, geometry or poni
- Returns:
same geometry as a Fit2dGeometry named-tuple
- pyFAI.geometry.fit2d.degrees(x)
- pyFAI.geometry.fit2d.radians(x)