condor.utils package¶
Submodules¶
condor.utils.bodies module¶
-
condor.utils.bodies.get_icosahedron_normal_vectors()[source]¶ Return array of normal vectors at the faces of a regular icosahedron
-
condor.utils.bodies.get_icosahedron_vertices()[source]¶ Return array of vertices vectors of a regular icosahedron
-
condor.utils.bodies.make_icosahedron_map(N, nRmax, extrinsic_rotation=None)[source]¶ Generate map of a uniform icosahedron (density = 1) on a regular grid
Orientation: The cartesian grid axis all lie parallel to 2-fold symmetry axes of the icosahedron.
- Args:
N (int): Edge length of the grid in unit pixels nRmax (float): Outer radius of the icosahedron in unit pixels - Kwargs:
- :rotation (
condor.utils.rotation.Rotation): Rotation instance for extrinsic rotation of the icosahedron.
-
condor.utils.bodies.make_icosahedron_map_slow(N, nRmax, extrinsic_rotation=None)[source]¶ Generate map of a uniform icosahedron (density = 1) on a regular grid (slow python implementation)
Orientation: The cartesian grid axis all lie parallel to 2-fold symmetry axes of the icosahedron.
- Args:
N (int): Edge length of the grid in unit pixels nRmax (float): Outer radius of the icosahedron in unit pixels
Kwargs:
:rotation (condor.utils.rotation.Rotation): Rotation instance for extrinsic rotation of the icosahedron.
-
condor.utils.bodies.make_sphere_map(N, nR)[source]¶ Generate a 3D map of a sphere particle on a regular grid (values between 0 and 1)
The result is quite rough (i.e. linear interpolation)
- Args:
N (int): Edge length of the grid in unit pixels nR (float): Radius in unit pixels
Note
This function was written for testing purposes and generates a map with rough edges. Use
condor.particle.particle_sphere.ParticleSpherefor more accurate uniform sphere diffraction simulations.
-
condor.utils.bodies.make_spheroid_map(N, nA, nC, rotation=None)[source]¶ Generate a 3D binary map of a spheroid particle on a regular grid
The result is very rough (i.e. nearest-neighbor interpolation)
- Args:
N (int): Edge length of the grid in unit pixels nA (float): Radius perpendicular to the rotation axis of the ellipsoid in unit pixels nC (float): Radius along the rotation axis of the ellipsoid in unit pixels
Kwargs:
:rotation (condor.utils.rotation.Rotation): Rotation instance for extrinsic rotation of the icosahedron.Note
This function was written for testing purposes and generates a map with rough edges. Use
condor.particle.particle_spheroid.ParticleSpheroidfor more accurate uniform spheroid diffraction simulations.
condor.utils.config module¶
Reading, writing and converting configuration in different representations
condor.utils.diffraction module¶
-
condor.utils.diffraction.crystallographic_resolution(wavelength, pixel_center_distance, detector_distance)[source]¶ Returns crystallographic resolution \(R_f\) (full-period resolution) in unit meter
\[R_f = \frac{ \lambda }{ 2 \sin\left( \arctan\left( \frac{X}{D} \right) / 2 \right) }\]- Args:
wavelength (float): Photon wavelength \(\lambda\) in unit meter pixel_center_distance (float): Distance \(X\) between beam center and pixel measured orthogonally with respect to the beam axis. The unit is meter detector_distance: Distance \(D\) between interaction point and detector plane in unit meter
-
condor.utils.diffraction.nyquist_pixel_size(wavelength, detector_distance, particle_size)[source]¶ Returns size \(p_N\) of one Nyquist pixel on the detector in unit meter
\[p_N = \frac{ D \lambda }{ d }\]- Args:
wavelength (float): Photon wavelength \(\lambda\) in unit meter detector_distance (float): Distance \(D\) between interaction point and detector plane in unit meter particle_size (float): Size or characteristic dimension \(d\) of the particle in unit meter
-
condor.utils.diffraction.polarization_factor(x, y, detector_distance, polarization='ignore')[source]¶ Returns polarization factor for a given geometry and polarization
Horizontally polarized:
\[P = \cos^2\left(rcsin\left(\]rac{x}{sqrt{x^2+y^2+D^2}} ight) ight)
Vertically polarized:
\[P = \cos^2\left(rcsin\left(\]rac{y}{sqrt{x^2+y^2+D^2}} ight) ight)
Unpolarized:
P = left(1 + cos^2left(rac{sqrt{x^2+y^2}}{sqrt{x^2+y^2+D^2}} ight)^2 ight)
Ignore polarization:
\[P = 1\]- Args:
x (float): horizontal pixel coordinate \(x\) with respect to beam center in unit meter y (float): vertical pixel coordinate \(y\) with respect to beam center in unit meter detector_distance (float): detector distance \(D\) in unit meter polarization (str): Type of polarization can be either vertical, horizontal, unpolarized, or ignore
-
condor.utils.diffraction.resolution_element(wavelength, pixel_center_distance, detector_distance)[source]¶ Returns length \(R_h\) of one resolution element (half-period resolution) in unit meter
\[R_h = \frac{ \lambda }{ 4 \, \sin\left( \arctan \left( \frac{X}{D} \right) / 2 \right) }\]- Args:
wavelength (float): Photon wavelength \(\lambda\) in unit meter pixel_center_distance (float): Distance \(X\) between beam center and pixel measured orthogonally with respect to the beam axis. The unit is meter detector_distance: Distance \(D\) between interaction point and detector plane in unit meter
condor.utils.linalg module¶
condor.utils.log module¶
-
condor.utils.log.log_and_raise_error(logger, message)¶
-
condor.utils.log.log_debug(logger, message)¶
-
condor.utils.log.log_info(logger, message)¶
-
condor.utils.log.log_warning(logger, message)¶
condor.utils.material module¶
-
class
condor.utils.material.AbstractMaterial[source]¶ -
get_attenuation_length(photon_wavelength)[source]¶ Return the absorption length in unit meter for the given wavelength \(\lambda\)
\[\mu = \frac{1}{\rho\,\mu_a(\lambda)}\]\(\rho\): Average atom density
\(\mu_a(\lambda)\): Photoabsorption cross section at photon energy \(\lambda\)
- Args:
photon_wavelength (float): Photon wavelength in unit meter
[Henke1993] B.L. Henke, E.M. Gullikson, and J.C. Davis. X-ray interactions: photoabsorption, scattering, transmission, and reflection at E=50-30000 eV, Z=1-92, Atomic Data and Nuclear Data Tables Vol. 54 (no.2), 181-342 (July 1993).
See also http://henke.lbl.gov/.
-
get_beta(photon_wavelength)[source]¶ Return \(\beta\) (imaginary part of \(\delta n\)) at a given wavelength
\[n = 1 - \delta n = 1 - \delta - i \beta\]\(n\): Refractive index
See also
condor.utils.material.Material.get_n()- Args:
photon_wavelength (float): Photon wavelength in unit meter
-
get_delta(photon_wavelength)[source]¶ Return \(\delta\) (real part of \(\delta n\)) at a given wavelength
\(n\): Refractive index
See also
condor.utils.material.Material.get_n()- Args:
photon_wavelength (float): Photon wavelength in unit meter
-
get_dn(photon_wavelength)[source]¶ Return \(\delta n\) at a given wavelength
\[\delta n = 1 - n\]\(n\): Refractive index
See also
condor.utils.material.Material.get_n()- Args:
photon_wavelength (float): Photon wavelength in unit meter
-
get_n(photon_wavelength)[source]¶ Return complex refractive index at a given wavelength (Henke, 1994)
\[n = 1 - \frac{ r_0 }{ 2\pi } \lambda^2 \sum_i \rho_i f_i(0)\]\(r_0\): classical electron radius
\(\rho_q\): atomic number density of atom species \(i\)
\(f_q(0)\): atomic scattering factor (forward scattering) of atom species \(i\)
- Args:
photon_wavelength (float): Photon wavelength in unit meter
-
get_photoabsorption_cross_section(photon_wavelength)[source]¶ Return the photoabsorption cross section \(\mu_a\) at a given wavelength \(\lambda\)
\[\mu_a = 2 r_0 \lambda f_2\]\(r_0\): classical electron radius
\(f_2\): imaginary part of the atomic scattering factor
- Args:
photon_wavelength (float): Photon wavelength in unit meter
-
get_transmission(thickness, photon_wavelength)[source]¶ Return transmission coefficient \(T\) for given material thickness \(t\) and wavelength \(\lambda\) [Henke1993]
\[T = e^{-\rho\,\mu_a(\lambda)\,t}\]\(\rho\): Average atom density
\(\mu_a(\lambda)\): Photoabsorption cross section at photon energy \(\lambda\)
- Args:
thickness (float): Material thickness in unit meter photon_wavelength (float): Photon wavelength in unit meter
[Henke1993] B.L. Henke, E.M. Gullikson, and J.C. Davis. X-ray interactions: photoabsorption, scattering, transmission, and reflection at E=50-30000 eV, Z=1-92, Atomic Data and Nuclear Data Tables Vol. 54 (no.2), 181-342 (July 1993).
See also http://henke.lbl.gov/.
-
-
class
condor.utils.material.AtomDensityMaterial(material_type, massdensity=None, atomic_composition=None)[source]¶ Bases:
condor.utils.material.AbstractMaterialClass for material model
- Args:
material_type (str): The material type can be either customor one of the standard types, i.e. tabulated combinations of massdensity and atomic composition, listed herecondor.utils.material.MaterialType.- Kwargs:
massdensity (float): Mass density in unit kilogram per cubic meter (default None)atomic_composition (dict): Dictionary of key-value pairs for atom species (e.g. 'H'for hydrogen) and concentration (defaultNone)
-
get_atomic_composition(normed=False)[source]¶ Return dictionary of atomic concentrations
- Args:
normed (bool): If Truethe concentrations are rescaled by a common factor such that their sum equals 1 (defaultFalse)
-
get_electron_density()[source]¶ Return electron density \(\rho_e\) in unit inverse cubic meter
\[\rho_e = \frac{\rho_m \cdot \sum_i}{\left( \sum_i c_i m_i \right) \left( \sum_i c_i Z_i \right)}\]\(\rho_m\): Mass denisty of material
\(c_i\): Normalised fraction of atom species \(i\)
\(m_i\): Standard atomic mass of atom species \(i\)
\(Z_i\): Atomic number of atom species \(i\)
-
get_f(photon_wavelength)[source]¶ Get effective average complex scattering factor for forward scattering at a given photon wavlength from Henke tables
- Args:
photon_wavlength (float): Photon wavelength in unit meter
-
get_scatterer_density()[source]¶ Return total atom density \(\rho\) in unit inverse cubic meter
\[\rho = \frac{\rho_m}{\sum_i c_i m_i}\]\(\rho_m\): Mass denisty of material
\(c_i\): Normalised fraction of atom species \(i\)
\(m_i\): Standard atomic mass of atom species \(i\)
-
set_atomic_concentration(element, relative_concentration)[source]¶ Set the concentration of a given atomic species
- Args:
element (str): Atomic species (e.g. 'H'for hydrogen)relative_concentration (float): Relative quantity of atoms of the given atomic species with respect to the others (e.g. for water: hydrogen concentration 2., oxygen concentration1.)
-
class
condor.utils.material.ElectronDensityMaterial(electron_density)[source]¶ Bases:
condor.utils.material.AbstractMaterialClass for electron density material model
Thomson scattering with the given value for the electron density is used to determine the material’s scattering properties.
- Args:
electron_density: (float): Electron density in unit inverse cubic meter
-
class
condor.utils.material.MaterialType[source]¶ Standard material types:
material_type\(\rho_m\) [kg/m3] Atomic composition Reference custommassdensityatomic_composition'water'995 (25 deg. C) \(H_2O\) [ONeil1868] p. 1868 'protein'1350 \(H_{86}C_{52}N_{13}O_{15}S\) [Bergh2008] 'dna'1700 \(H_{11}C_{10}N_4O_6P\) [Bergh2008] 'lipid'1000 \(H_{69}C_{36}O_6P\) [Bergh2008] 'cell'1000 \(H_{23}C_3NO_{10}S\) [Bergh2008] 'poliovirus'1340 \(C_{332652}H_{492388}N_{98245}O_{131196}P_{7501}S_{2340}\) [Molla1991] 'styrene'902 (25 deg. C) \(C_8H_8\) [Haynes2013] p. 3-488 'sucrose'1581 (17 deg. C) \(C_{12}H_{22O1}\) [Lide1998] p. 3-172 -
atomic_compositions= {'water': {'C': 0.0, 'H': 2.0, 'O': 1.0, 'N': 0.0, 'P': 0.0, 'S': 0.0}, 'cell': {'C': 3.0, 'H': 23.0, 'O': 10.0, 'N': 1.0, 'P': 0.0, 'S': 1.0}, 'dna': {'C': 10.0, 'H': 11.0, 'O': 6.0, 'N': 4.0, 'P': 1.0, 'S': 0.0}, 'lipid': {'C': 36.0, 'H': 69.0, 'O': 6.0, 'N': 0.0, 'P': 1.0, 'S': 0.0}, 'poliovirus': {'C': 332652.0, 'H': 492388.0, 'O': 131196.0, 'N': 98245.0, 'P': 7501.0, 'S': 2340.0}, 'styrene': {'C': 8.0, 'H': 8.0, 'O': 0.0, 'N': 0.0, 'P': 0.0, 'S': 0.0}, 'protein': {'C': 52.0, 'H': 86.0, 'O': 15.0, 'N': 13.0, 'P': 0.0, 'S': 1.0}, 'sucrose': {'C': 12.0, 'H': 22.0, 'O': 11.0, 'N': 0.0, 'P': 0.0, 'S': 0.0}}¶ Dictionary of atomic compositions (available keys are the tabulated
material_types)
-
mass_densities= {'water': 995.0, 'cell': 1000.0, 'dna': 1700.0, 'lipid': 1000.0, 'poliovirus': 1340.0, 'styrene': 902.0, 'protein': 1350.0, 'sucrose': 1581.0}¶ Dictionary of mass densities (available keys are the tabulated
material_types)
-
-
condor.utils.material.atomic_names= ['H', 'He', 'Li', 'Be', 'B', 'C', 'N', 'O', 'F', 'Ne', 'Na', 'Mg', 'Al', 'Si', 'P', 'S', 'Cl', 'Ar', 'K', 'Ca', 'Sc', 'Ti', 'V', 'Cr', 'Mn', 'Fe', 'Co', 'Ni', 'Cu', 'Zn', 'Ga', 'Ge', 'As', 'Se', 'Br', 'Kr', 'Rb', 'Sr', 'Y', 'Zr', 'Nb', 'Mo', 'Tc', 'Ru', 'Rh', 'Pd', 'Ag', 'Cd', 'In', 'Sn', 'Sb', 'Te', 'I', 'Xe', 'Cs', 'Ba', 'La', 'Ce', 'Pr', 'Nd', 'Pm', 'Sm', 'Eu', 'Gd', 'Tb', 'Dy', 'Ho', 'Er', 'Tm', 'Yb', 'Lu', 'Hf', 'Ta', 'W', 'Re', 'Os', 'Ir', 'Pt', 'Au', 'Hg', 'Tl', 'Pb', 'Bi', 'Po', 'At', 'Rn', 'Fr', 'Ra', 'Ac', 'Th', 'Pa', 'U', 'Np', 'Pu', 'Am', 'Cm', 'Bk', 'Cf', 'Es', 'Fm', 'Md', 'No', 'Lr', 'Rf', 'Db', 'Sg', 'Bh', 'Hs', 'Mt', 'Ds', 'Rg', 'Cp', 'Uut', 'Uuq', 'Uup', 'Uuh', 'Uus', 'Uuo']¶ List of atom names (i.e. latin abbreviations) for all elements sorted by atomic number (increasing order).
-
condor.utils.material.atomic_numbers= [1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23, 24, 25, 26, 27, 28, 29, 30, 31, 32, 33, 34, 35, 36, 37, 38, 39, 40, 41, 42, 43, 44, 45, 46, 47, 48, 49, 50, 51, 52, 53, 54, 55, 56, 57, 58, 59, 60, 61, 62, 63, 64, 65, 66, 67, 68, 69, 70, 71, 72, 73, 74, 75, 76, 77, 78, 79, 80, 81, 82, 83, 84, 85, 86, 87, 88, 89, 90, 91, 92, 93, 94, 95, 96, 97, 98, 99, 100, 101, 102, 103, 104, 105, 106, 107, 108, 109, 110, 111, 112, 113, 114, 115, 116, 117, 118]¶ List of atomic numbers of all elements in increasing order.
-
condor.utils.material.get_atomic_mass(element)¶ Returns the atomic mass (standard atomic weight in unit Dalton) for a given element.
- Args:
element (str): Element name (abbreviation of the latin name, for example ‘He’ for helium).
-
condor.utils.material.get_atomic_number(element)¶ Returns the atomic number for a given element.
- Args:
element (str): Element name (abbreviation of the latin name, for example ‘He’ for helium).
-
condor.utils.material.get_atomic_scattering_factors(element)¶ Returns 2-dim. array of photon energy [eV] vs. real and imaginary part of the atomic scattering factor (forward scattering) for a given element.
- Args:
element (str): Element name (abbreviation of the latin name, for example ‘He’ for helium).
condor.utils.photon module¶
-
class
condor.utils.photon.Photon(wavelength=None, energy=None, energy_eV=None, frequency=None)[source]¶ Class for X-ray photon
A Photon instance is initialised with one keyword argument - either with
wavelength,energyorenergy_eV.- Kwargs:
wavelength (float): Photon wavelength in unit meter energy (float): Photon energy in unit Joule energy_eV (float): Photon energy in unit electron volt frequency (float): Photon frequency in unit Hz
-
set_energy(energy)[source]¶ Set photon energy in unit Joule
- Args:
energy (float): Photon energy in unit Joule
-
set_energy_eV(energy_eV)[source]¶ Set photon energy in unit electron volt
- Args:
energy (float): Photon energy in unit electron volt
condor.utils.pixelmask module¶
-
class
condor.utils.pixelmask.PixelMask[source]¶ CXI bitmasks follow the formalism described in Maia (2012) [1] (see also www.cxidb.org) and encode the status of detector pixels. In Condor these bitmasks are represented by arrays of unsigned 16 bit integer values.
The status of a pixel is a combination of pixel properties represented by the bits of the bitmask value. A pixel with no bit set (all bits are 0) represents an “perfect” pixel. No property is set.
PIXEL_IS_PERFECT= 0Each bit that is set to 1 indicates that the respective property is assigned to that pixel.
Bit code:
Name Bit position Bit in integer representation PIXEL_IS_INVALID0th 0 PIXEL_IS_SATURATED1st 2 PIXEL_IS_HOT2nd 4 PIXEL_IS_DEAD3rd 8 PIXEL_IS_SHADOWED4th 16 PIXEL_IS_IN_PEAKMASK5th 32 PIXEL_IS_TO_BE_IGNORED6th 64 PIXEL_IS_BAD7th 128 PIXEL_IS_OUT_OF_RESOLUTION_LIMITS8th 256 PIXEL_IS_MISSING9th 512 PIXEL_IS_NOISY10th 1024 PIXEL_IS_ARTIFACT_CORRECTED11th 2048 PIXEL_FAILED_ARTIFACT_CORRECTION12th 4096 PIXEL_IS_PEAK_FOR_HITFINDER13th 8192 PIXEL_IS_PHOTON_BACKGROUND_CORRECTED14th 16384 (not used) 15th 32768 [1] Filipe R. N. C. Maia, The Coherent X-ray Imaging Data Bank, Nature Methods 9, 854-855 (2012) doi:10.1038/nmeth.2110. -
PIXEL_FAILED_ARTIFACT_CORRECTION= 4096¶
-
PIXEL_IS_ARTIFACT_CORRECTED= 2048¶
-
PIXEL_IS_BAD= 128¶
-
PIXEL_IS_DEAD= 8¶
-
PIXEL_IS_HOT= 4¶
-
PIXEL_IS_INVALID= 1¶
-
PIXEL_IS_IN_MASK= 767¶
-
PIXEL_IS_IN_PEAKMASK= 32¶
-
PIXEL_IS_MISSING= 512¶
-
PIXEL_IS_NOISY= 1024¶
-
PIXEL_IS_OUT_OF_RESOLUTION_LIMITS= 256¶
-
PIXEL_IS_PEAK_FOR_HITFINDER= 8192¶
-
PIXEL_IS_PERFECT= 0¶
-
PIXEL_IS_PHOTON_BACKGROUND_CORRECTED= 16384¶
-
PIXEL_IS_SATURATED= 2¶
-
PIXEL_IS_SHADOWED= 16¶
-
PIXEL_IS_TO_BE_IGNORED= 64¶
-
condor.utils.profile module¶
-
class
condor.utils.profile.Profile(model, focus_diameter)[source]¶ Class for spatial illumination profile
Arguments:
model (str): See condor.utils.profile.Profile.set_model()focus_diameter (float): Focus diameter / full-width half maximum in unit meter -
set_model(model)[source]¶ Set the model
Args:
model (str): Model for the spatial illumination profile
Choose one of the following options:
None- flat illumination profile of infinite extent and the intensity at every point corresponds to the intensity of a top hat profile (see below)'top_hat'- Circular top hat profile'pseudo_lorentzian'- 2D radially symmetrical profile composed of two Gaussian profiles obtained by a fit to a Lorentzian profile (used instead of a real Lorentzian because this function can be normalised)'gaussian'- 2D radially symmetrical Gaussian profile
-
condor.utils.resample module¶
-
condor.utils.resample.downsample(array2d0, factor0, mode='pick', mask2d0=None, bad_bits=None, min_N_pixels=1)[source]¶
-
condor.utils.resample.downsample_pos(pos, size, binning)¶
-
condor.utils.resample.upsample_pos(pos, size, binning)¶
condor.utils.rotation module¶
This module is an implementation of a variety of tools for rotations in 3D space.
-
condor.utils.rotation.R_x(t)¶
-
condor.utils.rotation.R_y(t)¶
-
condor.utils.rotation.R_z(t)¶
-
class
condor.utils.rotation.Rotation(values=None, formalism=None)[source]¶ Class for a rotation in 3D space
Arguments:
values (array): Array of values that define the rotation. For random rotations set values=``None`` and for example formalism=``’random’
. (default ``None)formalism: Formalism that defines how the argument values is interpreted. If
Noneno rotation. (defaultNone)Rotation formalism can be one of the following:
formalismVariables values'quaternion'\(q = w + ix + jy + kz\) \([w,x,y,z]\) 'rotation_matrix'\(R = \begin{pmatrix} R_{11} & R_{12} & R_{13} \\ R_{21} & R_{22} & R_{23} \\ R_{31} & R_{32} & R_{33} \end{pmatrix}\) \([[R_{11},R_{12},R_{13}],[R_{21},R_{22},R_{23}],[R_{31},R_{32},R_{33}]]\) 'euler_angles_zxz'\(e_1^{(z)}\), \(e_2^{(x)}\), \(e_3^{(z)}\) \([e_1^{(y)},e_2^{(z)},e_3^{(y)}]\) 'euler_angles_xyx'\(e_1^{(x)}\), \(e_2^{(y)}\), \(e_3^{(x)}\) \([e_1^{(x)},e_2^{(y)},e_3^{(x)}]\) 'euler_angles_xyz'\(e_1^{(x)}\), \(e_2^{(y)}\), \(e_3^{(z)}\) \([e_1^{(x)},e_2^{(y)},e_3^{(z)}]\) 'euler_angles_yzx'\(e_1^{(y)}\), \(e_2^{(z)}\), \(e_3^{(x)}\) \([e_1^{(y)},e_2^{(z)},e_3^{(x)}]\) 'euler_angles_zxy'\(e_1^{(z)}\), \(e_2^{(x)}\), \(e_3^{(y)}\) \([e_1^{(z)},e_2^{(x)},e_3^{(y)}]\) 'euler_angles_zyx'\(e_1^{(z)}\), \(e_2^{(y)}\), \(e_3^{(x)}\) \([e_1^{(z)},e_2^{(y)},e_3^{(x)}]\) 'euler_angles_yxz'\(e_1^{(y)}\), \(e_2^{(x)}\), \(e_3^{(z)}\) \([e_1^{(y)},e_2^{(x)},e_3^{(z)}]\) 'euler_angles_xzy'\(e_1^{(x)}\), \(e_2^{(z)}\), \(e_3^{(y)}\) \([e_1^{(x)},e_2^{(z)},e_3^{(y)}]\) 'random'fully random rotation None'random_x'random rotation around \(x\) axis None'random_y'random rotation around \(y\) axis None'random_z'random rotation around \(z\) axis None-
get_as_euler_angles(rotation_axes='zxz')[source]¶ Get rotation in Euler angle represantation \([e_1^{(z)}, e_2^{(x)}, e_3^{(z)}]\) (for the case of
rotation_axis='zxz').- Kwargs:
rotation_axes (str): Rotation axes of the three rotations (default 'zxz')
-
get_as_quaternion(unique_representation=False)[source]¶ Get rotation in quaternion representation \([w, x, y, z]\).
- Kwargs:
unique_representation (bool): Make quaternion unique. For more details see the documentation of condor.utils.rotation.unique_representation_quat()(default = False)
-
is_similar(rotation, tol=1e-05)[source]¶ Compare rotation with another instance of the Rotation class. If quaternion distance is smaller than tol return
True- Args:
- :rotation (
condor.utils.rotation.Rotation): Instance of the Rotation class - Kwargs:
tol (float): Tolerance for similarity. This is the maximum distance of the two quaternions in 4D space that will be interpreted for similar rotations. (default 0.00001)
-
rotate_vector(vector, order='xyz')[source]¶ Return the rotated copy of a given vector
- Args:
vector (array): 3D vector - Kwargs:
order (str): Order of geometrical axes in array representation of the given vector (default 'xyz')
-
rotate_vectors(vectors, order='xyz')[source]¶ Return the rotated copy of a given array of vectors
- Args:
vectors (array): Array of 3D vectors with shape (\(N\), 3) with \(N\) denoting the number of 3D vectors - Kwargs:
order (str): Order of geometrical axes in array representation of the given vector (default 'xyz')
-
set_with_euler_angles(euler_angles, rotation_axes='zxz')[source]¶ Set rotation with an array of three euler angles
- Args:
euler_angles: Array of the three euler angles representing consecutive rotations - Kwargs:
rotation_axes(str): Rotation axes of the three consecutive rotations (default = ‘zxz’)
-
set_with_quaternion(quaternion)[source]¶ Set rotation with a quaternion
- Args:
quaternion: Numpy array representing the quaternion \(w+ix+jy+kz\):
[\(w\), \(x\), \(y\), \(z\)] = [\(\cos(\theta/2)\), \(u_x \sin(\theta/2)\), \(u_y \sin(\theta/2)\), \(u_z \sin(\theta/2)\)]
with \(\theta\) being the rotation angle and \(\vec{u}=(u_x,u_y,u_z)\) the unit vector that defines the axis of rotation.
-
-
class
condor.utils.rotation.Rotations(values=None, formalism=None)[source]¶ Class for a list of rotations in 3D space
- Args:
values (array): Arrays of values that define the rotation. For random rotations set values = None(defaultNone)formalism (str): See condor.utils.rotation.Rotation. For no rotation setformalism = None(defaultNone)
-
condor.utils.rotation.euler_from_quat(q, rotation_axes='zxz')[source]¶ Return euler angles from quaternion (ref. [euclidianspace_mxToQuat], [euclidianspace_quatToEul]).
- Args:
q: Numpy array \([w,x,y,z]\) that represents the quaternion - Kwargs:
rotation_axes(str): Rotation axes of the three consecutive Euler rotations (default \'zxz\')
-
condor.utils.rotation.euler_to_rotmx(euler_angles, rotation_axes='zxz')[source]¶ Obtain rotation matrix from three euler angles and the rotation axes
- Args:
euler_angles (array): Length-3 array of euler angles - Kwargs:
rotation_axes (str): Rotation axes of the three consecutive Euler rotations (default 'zxz')
-
condor.utils.rotation.make_euler_unique_repax(euler)[source]¶ Make Euler angles for roations with a repeated axis (zxz, zyz, yzy, yxy, xzx, xyx) unique (such that three euler angles uniquely define a rotation).
There are the following ambiguities:
- Common circularity of rotations:
\[(e_1,e_2,e_3) \Leftrightarrow (e_1 + n \cdot 2\pi ,e_2 + n \cdot 2\pi ,e_3 + n \cdot 2\pi) \, , \, n \in \mathbb{Z}\]- Ambiguity resulting if first and last rotation is increased by \(\pi\):
\[(e_1,e_2,e_3) \Leftrightarrow (e_1 + \pi, -e_2, e_3 + \pi)\]- Gimbal lock if \((e_1+e_3) \mod 2\pi = 0\):
\[(e_1,e_2,e_3) \Leftrightarrow (a,e_2,-a)`, :math:`a \in \mathbb{Q}\]The follwing conventions are being used (observe order!):
- Allways:
\[(e_1 \, \textrm{mod} \, 2\pi, e_2 \, \textrm{mod} \, 2\pi, e_3 \, \textrm{mod} \, 2\pi)\]- If \(e_1 \geq \pi\):
\[((e_1 + pi) \mod 2\pi, 2\pi - e_2, (e_3 + \pi) \mod 2\pi)\]- If \(| e_1+e_3-2\pi | < \epsilon\):
\[(0,e_1,0)\]- Args:
- euler (array): Numpy array with the three euler angles in radian.
-
condor.utils.rotation.norm_quat(q, tolerance=1e-05)[source]¶ Return a copy of the normalised quaternion (adjust length to 1 if deviation larger than given tolerance).
- Args:
tolerance(float): Maximum deviation of length before rescaling (default 0.00001)
-
condor.utils.rotation.quat(theta, ux, uy, uz)¶
-
condor.utils.rotation.quat_conj(q)[source]¶ Return the conjugate quaternion as a length-4 array [w,-ix,-jy,-kz]
- Args:
q (array): Numpy array \([w,x,y,z]\) that represents the quaternion
-
condor.utils.rotation.quat_from_rotmx(R)[source]¶ Obtain the quaternion from a given rotation matrix (ref. [euclidianspace_mxToQuat])
- Args:
R: 3x3 array that represent the rotation matrix (see Conventions)
-
condor.utils.rotation.quat_mult(q1, q2)[source]¶ Return the product of two quaternions
- Args:
q1 (array): Length-4 array [\(w_1\), \(x_1\), \(y_1\), \(z_1\)] that represents the first quaternion q2 (array): Length-4 array [\(w_2\), \(x_2\), \(y_2\), \(z_2\)] that represents the second quaternion
-
condor.utils.rotation.quat_vec_mult(q, v)[source]¶ Return the product of a quaternion and a vector
- Args:
q (array): Length-4 array [\(w\), \(x\), \(y\), \(z\)] that represents the quaternion v (array): Length-3 array [\(v_x\), \(v_y\), \(v_z\)] that represents the vector
-
condor.utils.rotation.quat_x(theta)¶
-
condor.utils.rotation.quat_y(theta)¶
-
condor.utils.rotation.quat_z(theta)¶
-
condor.utils.rotation.rand_quat()[source]¶ Obtain a uniform random rotation in quaternion representation ([Shoemake1992] pages 129f)
-
condor.utils.rotation.rot_x(v, t)¶
-
condor.utils.rotation.rot_y(v, t)¶
-
condor.utils.rotation.rot_z(v, t)¶
-
condor.utils.rotation.rotate_quat(v, q)[source]¶ Return rotated version of a given vector by a given quaternion
- Args:
v (array): Length-3 array \([v_x,v_y,v_z]\) that represents the vector q (array): Length-4 array [\(w\), \(x\), \(y\), \(z\)] that represents the quaternion
-
condor.utils.rotation.rotmx_from_quat(q)[source]¶ Create a rotation matrix from given quaternion ([Shoemake1992] page 128)
- Args:
quaternion (array): \(q = w + ix + jy + kz\) ( values: \([w,x,y,z]\))
The direction of rotation follows the right hand rule
-
condor.utils.rotation.unique_representation_quat(q)[source]¶ Return the quaternion in a unique representation of rotations by avoiding the ambiguity that \(q\) and \(-q\) determine the same rotation.
The convention is that the first non-zero coordinate value of the quaternion array has to be positive.
- Args:
q (array): Length-4 array [\(w\), \(x\), \(y\), \(z\)] representing the qauternion \(q = w + ix + jy + kz\)
condor.utils.scattering_vector module¶
-
condor.utils.scattering_vector.generate_absqmap(X, Y, pixel_size, detector_distance, wavelength)[source]¶ Generate absolute scattering vector map from experimental parameters
- Args:
X (array): See condor.utils.scattering_vector.generate_qmap()Y (array): See condor.utils.scattering_vector.generate_qmap()pixel_size (float): See condor.utils.scattering_vector.generate_qmap()detector_distance (float): See condor.utils.scattering_vector.generate_qmap()wavelength (float): See condor.utils.scattering_vector.generate_qmap()
-
condor.utils.scattering_vector.generate_qmap(X, Y, pixel_size, detector_distance, wavelength, extrinsic_rotation=None, order='xyz')[source]¶ Generate scattering vector map from experimental parameters
- Args:
X (array): \(x\)-coordinates of pixels in unit meter Y (array): \(y\)-coordinates of pixels in unit meter pixel_size (float): Pixel size (i.e. edge length) in unit meter detector_distance (float): Distance from interaction point to detector plane wavelength (float): Photon wavelength in unit meter - Kwargs:
:extrinsic_rotation (
condor.utils.rotation.Rotation): Extrinsic rotation of the sample. IfNoneno rotation is applied (defaultNone)order (str): Order of scattering vector coordinates in the output array. Choose either 'xyz'or'zyx'(default'xyz')
-
condor.utils.scattering_vector.generate_qmap_3d(qn, qmax, extrinsic_rotation=None, order='xyz')[source]¶
-
condor.utils.scattering_vector.generate_rpix_3d(qn, qmax, wavelength, detector_distance, pixel_size)[source]¶
-
condor.utils.scattering_vector.q_from_p(p, wavelength)[source]¶ Return scattering vector from pixel position vector and photon wavelength.
- Args:
p (array): Pixel position vector \(\vec{p}=(p_x,p_y,p_z)\) with respect to the interaction point in unit meter wavelength (float): Photon wavelength in unit meter
condor.utils.sphere_diffraction module¶
-
condor.utils.sphere_diffraction.F_sphere_diffraction(K, q, r)¶ Scattering amplitude from homogeneous sphere (ref. [Feigin1987])
\[F(q) = \sqrt{K} \cdot f(q)\]\[f(q) = \frac{ 3 \left[ \sin(qr) - qr \cos(qr) \right]}{ (qr)^3 }\]\(I_0\): Primary intensity on the sample in unit number of photons per square meter
\(\rho_e\): Electron density in unit number of electrons per cubic meter
\(p\): Pixel size (i.e. edge length) in unit meter
\(D\): Detector distance in unit meter
\(r_0\): Classical electron radius in unit meter
\(V_r\): Sphere volume in unit cubic meter
- Args:
K (float): Intensity scaling factor \(K = I_0 \left(\rho_e \frac{p}{D} r_0 V_r\right)^2\) q (float/array): \(q\): Length of scattering vector in unit inverse meter r (float): \(r\): Sphere radius in unit meter
-
condor.utils.sphere_diffraction.I_sphere_diffraction(K, q, r)¶ Scattering Intensity from homogeneous sphere
\[I(q) = \left|F(q)\right|^2\]- Args:
K (float): See condor.utils.sphere_diffraction.F_sphere_diffraction()q (float/array): See condor.utils.sphere_diffraction.F_sphere_diffraction()r (float): \(r\): See condor.utils.sphere_diffraction.F_sphere_diffraction()
condor.utils.spheroid_diffraction module¶
-
condor.utils.spheroid_diffraction.F_spheroid_diffraction(K, q_x, q_y, a, c, theta, phi)¶ Scattering amplitude from homogeneous spheroid (ref. [Feigin1987], [Hamzeh1974])
The rotation axis is alligned parallel to the the \(y\)-axis before the rotations by \(theta\) and \(phi\) are applied (see below).
\[F(\vec{q}) = \sqrt{K} \cdot f(\vec{q})\]\[f(\vec{q}) = \frac{ 3 \left[\sin(qH(\vec{q})) - qH \cos(qH(\vec{q})) \right]}{ (qH(\vec{q}))^3}\]\[H(\vec{q}) = \sqrt{a^2 \sin^2(g(\vec{q}))+c^2 \cos^2(g(\vec{q}))}\]\[g(\vec{q}) = \arccos\left( \frac{ -q_x \cos(\theta) sin(\phi) + q_y \cos(\theta) cos(\phi) }{ q } \right)\]\(I_0\): Primary intensity on the sample in unit number of photons per square meter
\(\rho_e\): Electron density in unit number of electrons per cubic meter
\(p\): Pixel size (i.e. edge length) in unit meter
\(D\): Detector distance in unit meter
\(r_0\): Classical electron radius in unit meter
\(V_{a,c}\): Spheroid volume in unit cubic meter
- Args:
K (float): Intensity scaling factor \(K = I_0 \left(\rho_e \frac{p}{D} r_0 V_{a,c}\right)^2\) q_x (float/array): \(q_x\): \(x\)-coordinate of scattering vector in unit inverse meter q_y (float/array): \(q_y\): \(y\)-coordinate of scattering vector in unit inverse meter a (float): \(a\): radius perpendicular to the rotation axis of the ellipsoid in unit meter c (float): \(c\): radius along the rotation axis of the ellipsoid in unit meter theta (float): \(\theta\): extrinisc rotation around \(x\)-axis (1st, counter clockwise / right hand rule) phi (float): \(\phi\): extrinsic rotation around \(z\)-axis (2nd, counter clockwise / right hand rule)
-
condor.utils.spheroid_diffraction.I_spheroid_diffraction(K, qX, qY, a, c, theta, phi)¶ Scattering Intensity from homogeneous spheroid
\[I = \left|F\right|^2\]- Args:
K (float): See condor.utils.spheroid_diffraction.F_spheroid_diffraction()q_x (float/array): See condor.utils.spheroid_diffraction.F_spheroid_diffraction()q_y (float/array): See condor.utils.spheroid_diffraction.F_spheroid_diffraction()a (float): See condor.utils.spheroid_diffraction.F_spheroid_diffraction()c (float): See condor.utils.spheroid_diffraction.F_spheroid_diffraction()theta (float): See condor.utils.spheroid_diffraction.F_spheroid_diffraction()phi (float): See condor.utils.spheroid_diffraction.F_spheroid_diffraction()
-
condor.utils.spheroid_diffraction.to_spheroid_diameter(a, c)¶ Conversion from spheroid semi-diameters \(a\) and \(c\) to (sphere volume equivalent) diameter
-
condor.utils.spheroid_diffraction.to_spheroid_flattening(a, c)¶ Conversion from spheroid semi-diameters \(a\) and \(c\) to flattening (\(a/c\))
-
condor.utils.spheroid_diffraction.to_spheroid_semi_diameter_a(diameter, flattening)¶ Conversion from spheroid (sphere volume equivalent) diameter and flattening (\(a/c\)) to semi-diameter \(a\)
-
condor.utils.spheroid_diffraction.to_spheroid_semi_diameter_c(diameter, flattening)¶ Conversion from spheroid (sphere volume equivalent) diameter and flattening (\(a/c\)) to semi-diameter \(c\)
condor.utils.variation module¶
Tools for applying noise / statistical variation to values
-
class
condor.utils.variation.Variation(mode, spread, n=None, number_of_dimensions=1)[source]¶ Class for statistical variation of a variable
- Args:
mode (str): See condor.utils.variation.Variation.set_mode()spread (float/array): See condor.utils.variation.Variation.set_spread()- Kwargs:
n (int): Number of samples within the specified range (default None)number_of_dimensions (int): Number of dimensions of the variable (default 1)
-
get(v0)[source]¶ Get next value(s)
- Args:
v0 (float/int/array): Value(s) without variational deviation
-
get_conf()[source]¶ Get configuration in form of a dictionary. Another identically configured Variation instance can be initialised by:
conf = V0.get_conf() # V0: already existing Variation instance V1 = condor.utils.variation.Variation(**conf) # V1: new Variation instance with the same configuration as V0
-
set_mode(mode)[source]¶ Set the mode of variation
- Args:
mode (str): Mode of variation
Choose one of the following options
modeVariation model spreadparameterNoneNo variation 'uniform'Uniform random distribution Width 'normal'Normal random distribution Standard deviation 'poisson'Poissonian random distribution 'normal_poisson'Normal on top of Poissonian random dist. Std. dev. of normal distribution 'range'Equispaced grid around mean center position Width
-
set_spread(spread)[source]¶ Set spread of variation (standard deviation or full spread of values, see also
condor.utils.variation.Variation.set_mode())- Args:
spread (float): Width of the variation