reborn.simulate package¶
Submodules¶
reborn.simulate.clcore module¶
This module contains some core functions that are useful for simulating diffraction on GPU devices. Look to the reborn documentation to gain an understaning of how this module is meant to be used.
To get some information on compute devices (CPU/GPU) you can run the function clcore.help()
Some environment variables that affect the behaviour of this module:
REBORN_CL_GROUPSIZE : This sets the default groupsize.
PYOPENCL_CTX: This sets the device and platform automatically.
Using the above variables allows you to run the same code on different devices.
- reborn.simulate.clcore.get_all_gpu_devices()[source]¶
Search all platforms for GPU devices and return full list.
- Returns:
List of GPU devices. If none found, the list is empty.
- Return type:
- reborn.simulate.clcore.create_some_gpu_context()[source]¶
Since cl.create_some_context() sometimes forces a CPU on macs, this function will attempt to use a GPU context if possible. A CPU context will be returned if a GPU context is not found.
Returns: opencl context
- class reborn.simulate.clcore.ClCore(context=None, queue=None, group_size=32, double_precision=False, debug=0)[source]¶
Bases:
object
A container for the elementary building blocks that GPU diffraction simulations are composed of.
An instance of this class will initialize an opencl context and help maintain consistency in the device queue that compute jobs are sent to, along with consistency in the precision (double/single) when memory moves between CPU and GPU memory.
An instance of this class will attempt to help you manage an opencl context and command queue. You may choose the precision that you desire from the beginning, and this will be taken care of so that you don’t need to think about it when you move arrays between CPU and GPU memory.
An opencl context and queue will be created if you do not provide them. This is the most common mode of operation.
You may choose the group size, which is the number of compute units that have shared local memory. The environment variable REBORN_CL_GROUPSIZE can be used to set the default group size for a given machine if you want your code to be the same on multiple different machines.
The raw opencl kernel code will be compiled when an instance of this class is created.
- Parameters:
- get_device_name()[source]¶
Get the string (name) that describes the device that this ClCore instance is using. Example: ‘Intel(R) Gen9 HD Graphics NEO’
Note: You can find more properties of the device using e.g. ipython tab completion. This function simply does the following:
return self.context.devices[0].name
If you want a more complete list of device
Returns: string
- set_groupsize(group_size)[source]¶
If the environment variable REBORN_CL_GROUPSIZE is set then use that value.
If the group size exceeds the max allowed group size, then make it smaller (but print warning)
- vec4(x, dtype=None)[source]¶
Evdidently pyopencl does not deal with 3-vectors very well, so we use 4-vectors and pad with a zero at the end.
From Derek: I tested this at one point and found no difference… maybe newer pyopenCL is better..
This just does a trivial operation: return np.array([x.flat[0], x.flat[1], x.flat[2], 0.0], dtype=dtype)
- Parameters:
np.ndarray (x)
np.dtype (dtype) – Examples: np.complex, np.double
- Returns:
numpy array of length 4
- vec16(R, dtype=None)[source]¶
The best way to pass in a rotation matrix is as a float16. This is a helper function for preparing a numpy array so that it can be passed in as a float16.
From Derek: I had tested this and found no difference
See the vec4 function documentation also.
- Parameters:
R (numpy.ndarray) – input array
numpy.dtype (dtype) – default is np.float32
- Returns:
numpy array of length 16
- to_device(array=None, shape=None, dtype=None)[source]¶
This is a thin wrapper for pyopencl.array.to_device(). It will convert a numpy array into a pyopencl.array and send it to the device memory. So far this only deals with float and complex arrays, and it should figure out which type it is.
- Parameters:
array (numpy/cl array; float/complex type) – Input array.
shape (tuple) – Optionally specify the shape of the desired array. This is ignored if array is not None.
dtype (np.dtype) – Specify the desired type in opencl. The two types that are useful here are np.float32 and np.complex64
- Returns:
pyopencl array
- rotate_translate_vectors(rot, trans, vec_in, vec_out=None)[source]¶
Apply rotation followed by translation on GPU.
- Parameters:
rot – rotation matrix
trans – translation vector
v_in (Nx3 array) – input vectors
v_out (Nx3 array) – output vectors
Returns: output vectors
- test_rotate_vec(rot, trans, vec)[source]¶
Rotate a single vector. CPU arrays in, CPU array out. This is just for testing the consistency of memory allocation.
- mod_squared_complex_to_real(A, I, add=False)[source]¶
Compute the real-valued modulus square of complex numbers. Good example of a function that shouldn’t exist, but I needed to add it here because the pyopencl.array.Array class fails to do this operation correctly on some computers.
- Parameters:
A (clarray) – The complex amplitudes
I (clarray) – The real intensities
add (bool) – If false, overwrite the I array, else add to I
- Returns:
None
- sphere_form_factor(r, q, dens=None, a=None, add=False)[source]¶
Form factor \(f(q)\) for a sphere of radius \(r\), at given \(q\) magnitudes. The formula is
(1)¶\[f(q) = 4 \pi \frac{\sin(qr) - qr \cos(qr)}{q^3}\]When \(q = 0\), the following limit is used:
(2)¶\[f(0) = \frac{4}{3} \pi r^3\]Formula can be found, for example, in Table A.1 of Guinier (1964). There are no approximations in this formula beyond the 1st Born approximation; it is not a small-angle formula.
Note that you need to multiply this by the electron density of the sphere if you want reasonable amplitudes. E.g., water molecules have 10 electrons, a molecular weight of 18 g/mol and a density of 1 g/ml, so you can google search the electron density of water, which is 10*(1 g/cm^3)/(18 g/6.022e23) = 3.346e29 per m^3 .
- Parameters:
Returns: None
- phase_factor_qrf(q, r, f=None, R=None, U=None, a=None, add=False, twopi=False, n_chunks=1)[source]¶
Calculate diffraction amplitudes according to the sum:
(3)¶\[a_i = \sum_n f_n \exp(-i \vec{q}_i \cdot (\mathbf{R} \vec{r} + \vec{U}))\]- Parameters:
q (numpy/cl float array [N,3]) – Scattering vectors.
r (numpy/cl float array [M,3]) – Atomic coordinates.
f (numpy/cl complex array [M]) – Complex scattering factors.
R (numpy array [3,3]) – Rotation matrix.
U (numpy array) – Translation vector.
a (cl complex array [N]) – Complex scattering amplitudes (if you wish to manage your own opencl array).
add (bool) – True means add to the input amplitudes a rather than overwrite the amplitudes.
twopi (bool) – True means to multiply \(\vec{q}\) by \(2\pi\) just before calculating \(A(q)\).
n_r_chunks (int) – Run in n batches of position vectors to avoid memory issues.
- Returns:
Diffraction amplitudes. Will be a cl array if there are input cl arrays.
- Return type:
(numpy/cl complex array [N])
- phase_factor_pad(r, f=None, T=None, F=None, S=None, B=None, nF=None, nS=None, w=None, R=None, U=None, a=None, add=False, beam=None, pad=None)[source]¶
Calculate diffraction amplitudes according to the sum:
(4)¶\[a_i = \sum_n f_n \exp(-i \vec{q}_i \cdot (\mathbf{R} \vec{r} + \vec{U}))\]The \(\vec{q}_i\) vectors are computed on the GPU.
- Parameters:
r (numpy/cl float array [M,3]) – Atomic coordinates.
f (numpy/cl complex array [M]) – Complex scattering factors.
T – A 1x3 numpy array with vector components pointing from sample to the center of the first pixel in memory
F – A 1x3 numpy array containing the basis vector components pointing in the direction corresponding to contiguous pixels in memory (“fast scan”).
S – A 1x3 numpy array containing the basis vector components pointing in the direction corresponding to non-contiguous pixels in memory (“slow scan”).
B – A 1x3 numpy array with unit-vector components corresponding to the incident x-ray beam direction
nF – Number of fast-scan pixels (corresponding to F vector) in the detector panel
nS – Number of slow-scan pixels (corresponding to S vector) in the detector panel
w – The photon wavelength in meters
R – Optional numpy array [3x3] specifying rotation of atom vectors (we quietly transpose R and let it operate on q-vectors for speedups)
a – Optional output complex scattering amplitude cl array
- Returns:
A numpy array of length nF*nS containing complex scattering amplitudes
- Return type:
A
- phase_factor_mesh(r, f=None, N=None, q_min=None, q_max=None, dq=None, a=None, R=None, U=None, add=False, density_map=None, twopi=False)[source]¶
Compute the following sum on a regular 3D mesh of q samples:
(5)¶\[a_i = \sum_n f_n \exp(-i \vec{q}_i \cdot (\mathbf{R} \vec{r} + \vec{U}))\]The mesh is defined by the shape of the 3D array along with the minimum and maximum values of \(q_i\) along each edge \(i=1,2,3\). The vector components of q are computed according to: .. math:: q_{ni} = n Delta q_i + q_{text{min},i} where \(n\) is the array index (starting with \(n=0\)) along axis \(i\) and we define .. math:: Delta q_i = (q_{text{max},i} - q_{text{min},i})/(N_i - 1) for an array with shape \(N_i\).
- Parameters:
r (Nx3 numpy array) – Atomic coordinates
f (numpy array) – A numpy array of complex atomic scattering factors
N (numpy array length 3) – Number of q-space samples in each of the three dimensions
q_min (numpy array length 3) – Minimum q-space magnitudes in the 3d mesh. These values specify the center of the first voxel.
q_max (numpy array length 3) – Naximum q-space magnitudes in the 3d mesh. These values specify the center of the voxel.
dq (numpy array of length 3) – For legacy reasons, you can specify dq instead of N. Note that the relation between dq and N is dq = (q_max-q_min)/(N-1).
a (clArray) – device buffer, if available
R (3x3 array) – Rotation matrix, to be applied to r vectors
add (bool) – If True, add to device buffer a, else overwrite the buffer.
density_map (an instance of)
- Returns:
- An array of complex scattering amplitudes. By default this is a normal
numpy array. Optionally, this may be an opencl buffer.
- mesh_interpolation(a_map, q, N=None, q_min=None, q_max=None, dq=None, R=None, U=None, a=None, density_map=None, add=False, twopi=False)[source]¶
This is used in conjunction with
phase_factor_mesh
to interpolate amplitudes from a 3d mesh of simulated amplitudes.- Parameters:
a_map (numpy array) – Complex scattering amplitudes (usually generated from the function phase_factor_mesh()). Can also be real amplitudes.
N (int) – As defined in phase_factor_mesh()
q_min (float) – As defined in phase_factor_mesh()
q_max (float) – As defined in phase_factor_mesh()
q (Nx3 numpy array) – q-space coordinates at which we want to interpolate the complex amplitudes in a_dev
R (3x3 numpy array) – Rotation matrix that will act on the atom vectors (we quietly transpose R and let it operate on q-vectors for speedups)
a – (clarray) The output array (optional)
- Returns:
numpy array of complex amplitudes
- lattice_transform_intensities_pad(abc, N, T, F, S, B, nF, nS, w, R=None, I=None, add=False)[source]¶
Calculate crystal lattice transform intensities for a pixel-array detector. This is the usual transform for an idealized parallelepiped crystal (usually not very realistic…).
- Parameters:
abc (numpy array) – A 3x3 array containing real-space basis vectors. Vectors are contiguous in memory.
N (numpy array) – An array containing number of unit cells along each of three axes.
T (numpy array) – Translation to center of corner pixel.
F (numpy array) – Fast-scan basis vector.
S (numpy array) – Slow-scan basis vector.
B (numpy array) – Incident beam vector.
nF (int) – Number of fast-scan pixels.
nS (int) – Number of slow-scan pixels.
w (float) – Wavelength.
R (numpy array) – Rotation matrix acting on atom vectors. (we quietly transpose R and let it operate on q-vectors for speedups)
( (I) – class:pyopencl.array.Array) : OpenCL device array containing intensities.
add (bool) – If true, the function will add to the input I buffer, else the buffer is overwritten.
- Returns:
If I == None, then the output is a numpy array. Otherwise, it is an opencl array.
- gaussian_lattice_transform_intensities_pad(abc, N, T, F, S, B, nF, nS, w, R=None, I=None, add=False)[source]¶
Calculate crystal lattice transform intensities for a pixel-array detector. Uses a Gaussian approximation to the lattice transform.
- Parameters:
abc (numpy array) – A 3x3 array containing real-space basis vectors. Vectors are contiguous in memory.
N (numpy array) – An array containing number of unit cells along each of three axes.
T (numpy array) – Translation to center of corner pixel.
F (numpy array) – Fast-scan basis vector.
S (numpy array) – Slow-scan basis vector.
B (numpy array) – Incident beam vector.
nF (int) – Number of fast-scan pixels.
nS (int) – Number of slow-scan pixels.
w (float) – Wavelength.
R (numpy array) – Rotation matrix acting on atom vectors. (we quietly transpose R and let it operate on q-vectors for speedups)
( (I) – class:pyopencl.array.Array) : OpenCL device array containing intensities.
add (bool) – If true, the function will add to the input I buffer, else the buffer is overwritten.
- Returns:
If I == None, then the output is a numpy array. Otherwise, it is an opencl array.
- gaussian_lattice_transform_intensities(q, abc, N, R=None, I=None, add=False)[source]¶
Calculate crystal lattice transform intensities. Uses a Gaussian approximation to the lattice transform.
- Parameters:
q (
ndarray
) – q vectors.abc (
ndarray
) – A 3x3 array containing real-space basis vectors. Vectors are contiguous in memory.N (
ndarray
) – An array containing number of unit cells along each of three axes.R (
ndarray
) – Rotation matrix acting on atom vectors. (we quietly transpose R and let it operate on q-vectors for speedups)( (I) – class:pyopencl.array.Array) : OpenCL device array containing intensities.
add (bool) – If true, the function will add to the input I buffer, else the buffer is overwritten.
- Returns:
If I == None, then the output is a numpy array. Otherwise, it is an opencl array.
- reborn.simulate.clcore.phase_factor_qrf(q, r, f)[source]¶
Same as ClCore.phase_factor_qrf, but this skips the step of creating the ClCore instance first for convenience.
reborn.simulate.engines module¶
- class reborn.simulate.engines.MoleculePatternSimulator(beam=None, molecule=None, pad=None, oversample=10, max_mesh_size=200, clcore=None)[source]¶
Bases:
object
- class reborn.simulate.engines.SAXSPatternSimulator(pdb_file_path: str, beam: Beam, pad_geometry: PADGeometryList, concentration: float, thickness: float)[source]¶
Bases:
object
- Parameters:
pdb_file_path – Path to a PDB file.
beam – Beam information.
pad_geometry – Geometry information.
concentration – Protein concentration in kg/m^3.
thickness – Sample thickness in m.
- class reborn.simulate.engines.SAXSPatternFrameGetter(*args, **kwargs)[source]¶
Bases:
FrameGetter
When overriding the __init__ method, you should be sure to set the self.n_frames property. Do not use any positional or keyword arguments that cannot be serialized, such as file handlers (instead, provide a file path and create the handler within the __init__ method). If you pass in objects that cannot be serialized, your subclass will still work, but it will not be possible to pass your subclass into any of the parallelized processors that utilize FrameGetters
reborn.simulate.examples module¶
Some simple examples for testing purposes. Don’t build any of this into your code.
- reborn.simulate.examples.simulate_water(pad_geometry=None, beam=None, water_thickness=1e-06)[source]¶
Simulate water scatter. Takes a PAD geometry and beam specification and returns list of 2D numpy arrays with the scattering intensity in photon units.
- Parameters:
pad_geometry (list of
PADGeometry
’s) – List of PAD geometry specifications.beam (
Beam
) – Beam specificationwater_thickness (float) – Thickness of water in SI units.
- Returns:
list of 2D numpy arrays
- class reborn.simulate.examples.WaterFrameGetter(*args, **kwargs)[source]¶
Bases:
FrameGetter
Creates a
FrameGetter
based on the functionsimulate_water reborn.simulate.examples.simulate_water
. See the function arguments for class initialization (same arguments are accepted to initialize the class).When overriding the __init__ method, you should be sure to set the self.n_frames property. Do not use any positional or keyword arguments that cannot be serialized, such as file handlers (instead, provide a file path and create the handler within the __init__ method). If you pass in objects that cannot be serialized, your subclass will still work, but it will not be possible to pass your subclass into any of the parallelized processors that utilize FrameGetters
- thickness = None¶
- geom = None¶
- beam = None¶
- reborn.simulate.examples.lysozyme_molecule(pad_geometry=None, wavelength=1.5e-10, random_rotation=True)[source]¶
Simple simulation of lysozyme molecule using
ClCore
.- Parameters:
pad_geometry – List of
PADGeometry
instances.wavelength – As always, in SI units.
random_rotation – If True generate a random rotation. Default is True.
Returns: dictionary with {‘pad_geometry’: pads, ‘intensity’: data_list}
- class reborn.simulate.examples.PDBMoleculeSimulator(pdb_file=None, pad_geometry=None, beam=None, wavelength=1.5e-10, random_rotation=True)[source]¶
Bases:
object
A simple generator of simulated single-molecule intensities, from a pdb file. Defaults to lysozyme at 1.5 A wavelength, on a pnccd detector layout.
This will setup the opencl simulation core.
- Parameters:
pdb_file – path to a pdb file
pad_geometry – array of
PADGeometry
intanceswavelength – in SI units of course
random_rotation (bool) – True or False
- class reborn.simulate.examples.MoleculeSimulatorV1(beam=None, molecule=None, pad=None, oversample=10, max_mesh_size=200, clcore=None)[source]¶
Bases:
object
- class reborn.simulate.examples.CrystalSimulatorV1(pad_geometry=None, beam=None, crystal_structure=None, n_iterations=1, approximate_shape_transform=True, expand_symmetry=False, cl_double_precision=False, cl_group_size=32, poisson_noise=False)[source]¶
Bases:
object
Class for generating crystal diffraction patterns. Generates the average pattern upon:
Randomizing the x-ray beam direction (beam divergence)
Randomizing the outgoing beam according to pixel area (“pixel solid angle”)
Randomizing photon energy (spectral width)
Randomizing the orientation of crystal mosaic domains (mosaicity)
Randomizing the shape transforms or Gaussian crystal-size broadening
Computations are done on a GPU with OpenCL.
- Parameters:
pad_geometry (list of
PADGeometry
instances) – PAD geometry.beam (
Beam
) – A beam instance.crystal_structure (
CrystalStructure
) – A crystal structure.n_iterations (int) – Number of iterations to average over
approximate_shape_transform (bool) – Use a Gaussian approximation to shape transforms, else use the analytic parallelepiped shape transform formula.
expand_symmetry (bool) – Duplicate the asymmetric unit according to spacegroup symmetry in crystal_structure.
cl_double_precision (bool) – Use double precision if available on GPU device.
cl_group_size (int) – GPU group size (see the
ClCore
class).poisson_noise (bool) – Add Poisson noise to the resulting pattern.
- class reborn.simulate.examples.LysozymeFrameGetter(*args, **kwargs)[source]¶
Bases:
FrameGetter
When overriding the __init__ method, you should be sure to set the self.n_frames property. Do not use any positional or keyword arguments that cannot be serialized, such as file handlers (instead, provide a file path and create the handler within the __init__ method). If you pass in objects that cannot be serialized, your subclass will still work, but it will not be possible to pass your subclass into any of the parallelized processors that utilize FrameGetters
reborn.simulate.form_factors module¶
- reborn.simulate.form_factors.sphere_form_factor(radius, q_mags, check_divide_by_zero=True)[source]¶
Form factor \(f(q)\) for a sphere of radius \(r\), at given \(q\) magnitudes. The formula is
(6)¶\[f(q) = 4 \pi \frac{\sin(qr) - qr \cos(qr)}{q^3}\]When \(q = 0\), the following limit is used:
(7)¶\[f(0) = \frac{4}{3} \pi r^3\]Formula can be found, for example, in Table A.1 of Guinier (1964). There are no approximations in this formula beyond the 1st Born approximation; it is not a small-angle formula.
Note that you need to multiply this by the electron density of the sphere if you want reasonable amplitudes. E.g., water molecules have 10 electrons, a molecular weight of 18 g/mol and a density of 1 g/ml, so you can google search the electron density of water, which is 10*(1 g/cm^3)/(18 g/6.022e23) = 3.346e29 per m^3 .
- Parameters:
Returns: numpy array
- reborn.simulate.form_factors.spheroid_form_factor(R: ndarray, q_vecs: ndarray, check_divide_by_zero: bool = True)[source]¶
Form factor \(f(q)\) for a spheroid of radii \(a\), \(b\), \(c\), at given \(\vec{q}\) vectors. Assuming the vectors describing the semi-major axis vectors \(\vec{a}\), \(\vec{c}\), \(\vec{c}\) are contained in the rows of the matrix \(\text{R}\).
The form factor is
(8)¶\[\begin{split}f(\vec{q}) &= 4 \pi \text{T}(\mathrm{R}) \frac{\sin(h) - h \cos(h)}{h^3} \\ h &= | \mathrm{R}\vec{q} |^3\end{split}\]where \(\text{T}(\mathrm{R})\) is the product \(abc\) (i.e. the product of the eigenvalues of \(\text{T}(\mathrm{R})\)). When \(q = 0\), the following limit is used:
(9)¶\[f(0) = \frac{4}{3} \pi T(\mathrm{R})\]This formula may be derived by applying the Fourier stretch theorem to the case of a sphere. There are no approximations in this formula beyond the 1st Born approximation; it is not a small-angle formula.
Note that you need to multiply this by the electron density of the sphere if you want reasonable amplitudes. E.g., water molecules have 10 electrons, a molecular weight of 18 g/mol and a density of 1 g/ml, so you can google search the electron density of water, which is 10*(1 g/cm^3)/(18 g/6.022e23) = 3.346e29 per m^3 .
This code was generalized from code by Andrew Morgan, who credited Filipe Maia.
- Parameters:
R (np.ndarray) – Semi-major axes
q_vecs (numpy array) – Q vectors
check_divide_by_zero (bool) – Check for divide by zero. True by default.
Returns: numpy array
- reborn.simulate.form_factors.ellipsoid_form_factor(radius: float, eccentricity: float, theta: float, q_mags, check_divide_by_zero=True)[source]¶
Form factor \(f(q)\) for an ellipsoid of revolution (a ‘squashed’ or stretched’ sphere), at given \(q\) magnitudes. The formula is
(10)¶\[f(q) = 4\pi\epsilon \dfrac{\sin(qr_{\theta}) - qr_{\theta} \cos(qr_{\theta})}{q^3}\]Where
(11)¶\[r_{\theta} = r \sqrt{1 + (\epsilon^2 - 1) \cos^2{\theta}} = r \sqrt(\sin^2{\theta} + \epsilon^2 \cos^2{\theta})\]and \(theta\) is the angle with respect to the z-axis. When \(q = 0\), the following limit is used:
(12)¶\[f(0) = \frac{4}{3} \pi \epsilon r^3\]Formula can be found here: http://gisaxs.com/index.php/Form_Factor
Note that you need to multiply this by the electron density of the sphere if you want reasonable amplitudes. E.g., water molecules have 10 electrons, a molecular weight of 18 g/mol and a density of 1 g/ml, so you can google search the electron density of water, which is 10*(1 g/cm^3)/(18 g/6.022e23) = 3.346e29 per m^3 .
- Parameters:
radius (float) – In SI units of course.
eccentricity (float) – 1 for a sphere, less than 1 for an oblate spheroid and greater than 1 for a prolate spheroid.
theta (float) – angle with respect to the z-axis. Gives orientation of the ellipsoid.
q_mags (numpy array) – Also in SI units.
check_divide_by_zero (bool) – Check for divide by zero. True by default.
Returns: numpy array
reborn.simulate.gas module¶
- reborn.simulate.gas.isotropic_gas_intensity_profile(r_vecs=None, q_mags=None, atomic_numbers=None, photon_energy=None, molecule=None, beam=None)[source]¶
Calculate the isotropic scatter from one gas molecule using the Debye formula.
- Parameters:
r_vecs (
ndarray
) – Position vectors.q_mags (
ndarray
) – Q vectors.atomic_numbers (
ndarray
) – Atomic numbers.photon_energy (float) – Photon energy.
molecule (reborn.target.molecule.Molecule or str) – Molecule (overrides r_vecs and atomic_numbers)
beam (
Beam
) – Beam instance. Overrides photon_energy.
- Returns:
Intensity profile I(q)
- Return type:
- reborn.simulate.gas.get_gas_background(pad_geometry: PADGeometryList, beam: Beam, path_length=[0.0, 1.0], gas_type=None, temperature=293.15, pressure=101325.0, n_simulation_steps=20, poisson=False, verbose=False, as_list=False, __gamma=10, __polarization_off=False, __fixed_f2=None)[source]¶
Given a
PADGeometryList
along with aBeam
, calculate the scattering intensity \(I(q)\) of chamber gas.- Parameters:
pad_geometry (
PADGeometryList
) – PAD geometry info.beam (
Beam
) – X-ray beam parameters.path_length (list of float) – Start and end points corresponding to the gas that the beam passes through
gas_type (str, list) – String indicating the gas type, or list of strings indicating types
temperature (float) – Temperature of the gas.
pressure (float, list) – Pressure of the gas, or list of partial pressures
poisson (bool) – If True, add Poisson noise (default=True)
- Returns:
List of
ndarray
instances containing intensities.
reborn.simulate.solutions module¶
- reborn.simulate.solutions.get_water_profile(q, temperature=298)[source]¶
Depreciated: Use
water_scattering_factor_squared
- reborn.simulate.solutions.water_scattering_factor_squared(q, temperature=298, volume=None)[source]¶
Get water scattering profile \(|F(q)|^2\) via interpolations/smoothing of the data from Hura et al. 2003 (Physical Chemistry Chemical Physics 5(10), 2003,pp 1981-1991, https://doi.org/10.1039/B301481A) and Clark et al. 2010 (PNAS 107, 14003–14007). Temperatures from 1-77 C are included. Note that the returned values are the scattering factors \(|F(q)|^2\) per water molecule, so you should multiply this by the number of water molecules exposed to the x-ray beam. The number density of water is approximately 33.3679e27 / m^3. You may optionally provide the volume and the result will be multiplied by the number of water molecules.
- Parameters:
- Returns:
The scattering factor \(|F(q)|^2\) , possibly multiplied by the number of water molecules if you specified a volume.
- Return type:
- reborn.simulate.solutions.get_pad_solution_intensity(pad_geometry, beam, thickness=1e-05, liquid='water', temperature=298.0, poisson=True, as_list=True)[source]¶
Given a list of
PADGeometry
instances along with aBeam
instance, calculate the scattering intensity \(I(q)\) of a liquid of given thickness.note:
This function is only for convenience. Consider using :func:`water_scattering_factor_squared <reborn.simulate.solutions.water_scattering_factor_squared>` if speed is important; this will avoid re-calculating q magnitudes, solid angles, and polarization factors.
- Parameters:
pad_geometry (list of
PADGeometry
instances) – PAD geometry info.beam (
Beam
) – X-ray beam parameters.thickness (float) – Thickness of the liquid (assumed to be a sheet geometry)
liquid (str) – We can only do “water” at this time.
temperature (float) – Temperature of the liquid.
poisson (bool) – If True, add Poisson noise (default=True)
as_list (bool) – If True, return a list of 2D PAD arrays. Else return a flat 1D array. Default: True.
- Returns:
List of
ndarray
instances containing intensities.