reborn.target package

Submodules

reborn.target.atoms module

Resources for working with atoms: atomic scattering factors, electron densities.

reborn.target.atoms.atomic_symbols_to_numbers(symbols)[source]

Convert atomic symbol strings (such as C, N, He, Be) to atomic numbers.

Parameters:

symbols (string/list-like) – Atomic symbol strings. Case insensitive.

Returns:

Atomic numbers.

Return type:

numbers (ndarray)

reborn.target.atoms.atomic_numbers_to_symbols(numbers)[source]

Convert atomic numbers to atomic symbols (such as C, N, He, Ca)

Parameters:

numbers (ndarray) – Atomic numbers

Returns:

Atomic element symbols (such as H, He, Be)

Return type:

ndarray

reborn.target.atoms.henke_scattering_factors(atomic_numbers, photon_energies)[source]

Get complex atomic scattering factors from Henke tables.

Note: The inputs are converted to numpy arrays, and the output is a squeezed 2D array of shape (atomic_numbers.size, photon_energies.size).

Parameters:
  • atomic_numbers (int/list-like) – Atomic numbers.

  • photon_energies (float/list-like) – Photon energy in SI units.

Returns:

Complex scattering factors.

Return type:

ndarray

reborn.target.atoms.henke_dispersion_corrections(atomic_numbers, photon_energies)[source]

Same as henke_scattering_factors but subtracts the atomic number :param atomic_numbers: :param photon_energies:

Returns:

reborn.target.atoms.get_scattering_factors(atomic_numbers, photon_energy=None, beam=None)[source]

Get complex atomic scattering factors (from Henke tables) for a single photon energy and a range of atomic numbers. The scattering factor data come from the function get_henke_data. See the function get_scattering_factors_fixed_z if you have a range of photon energies and one atomic number.

Parameters:
  • atomic_numbers (int/list-like) – Atomic numbers.

  • photon_energy (float) – Photon energy in SI units.

  • beam (Beam) – Supply a beam instead of photon energy.

Returns:

Complex scattering factors.

Return type:

ndarray

reborn.target.atoms.get_scattering_factors_fixed_z(atomic_number, photon_energies)[source]

Get complex atomic scattering factors (from Henke tables) for a single atomic number and a range of photon energies.

Parameters:
  • atomic_number (int) – Atomic number.

  • photon_energies (float/list-like) – Photon energies in SI units.

Returns:

Complex scattering factors.

Return type:

ndarray

reborn.target.atoms.xraylib_scattering_factors(q_mags, atomic_number, photon_energy=None, beam=None)[source]

Get the q-dependent atomic scattering factors from the xraylib package. The scattering factor is equal to

\(\text{FF_Rayl}(Z, q) + \text{Fi}(Z, E) - i \text{Fii}(Z, E)\)

where FF_Rayl, Fi, and Fii are functions in the xraylib package. Note the strange minus sign in this formula, which is not a typo – it has been difficult to find documentation on these functions. Note that the dispersion relations in the Henke tables appear to be better than those in xraylib. See hubbel_henke_atomic_scattering_factors for a better alternative. Note that this function can be slow because the loops over q magnitudes are on the Python level because of the Python wrapper included with xraylib.

Parameters:
  • q_mags (ndarray) – q-vector magnitudes, where \(\vec{q} = \frac{2\pi}{\lambda}(\vec{s}-\vec{s}_0)\)

  • atomic_number (int) – The atomic number.

  • photon_energy (float) – The photon energy in SI units.

  • beam (Beam) – Supply a beam instead of photon energy.

Returns:

Complex scattering factors \(f(q)\)

Return type:

ndarray

reborn.target.atoms.xraylib_refractive_index(compound='H2O', density=1000, photon_energy=1.6021766339999998e-15, beam=None, approximate=False)[source]

Get the x-ray refractive index given a chemical compound, density, and photon energy.

Parameters:
  • compound (str) – Chemical compound formula (e.g. H20)

  • density (float) – Mass density in SI

  • photon_energy (float) – Photon energy in SI

  • beam (Beam) – For convenience, an alternative to photon_energy

  • approximate (bool) – Approximate with the non-resonant, high-frequency limit (Equation 1.1 of Guinier (1964)) (Default: False)

  • beam – Supply a beam instead of photon energy.

Returns:

Refractive index

Return type:

float

reborn.target.atoms.xraylib_scattering_density(compound='H2O', density=1000, photon_energy=1.6021766339999998e-15, beam=None, approximate=False)[source]

Calculate scattering density of compound given density and photon energy. This is approximately equal to the electron density far from resonance.

Parameters:
  • compound (str) – Chemical compound formula (e.g. H20)

  • density (float) – Mass density in SI

  • photon_energy (float) – Photon energy in SI

  • approximate (bool) – Approximate with the non-resonant, high-frequency limit (Equation 1.1 of Guinier (1964)) (Default: False)

  • beam (Beam) – Supply a beam instead of photon energy.

Returns:

Scattering density

Return type:

float

reborn.target.atoms.hubbel_form_factors(q_mags, atomic_number)[source]

Get the q-dependent atomic form factors. This allows for an arbitrary list of q magnitudes and returns an array. The scattering factors come from Hubbel et al 1975, and are accessed through the xraylib package.

Parameters:
  • q_mags (ndarray) – q vector magnitudes.

  • atomic_number (int) – Atomic number.

Returns:

Atomic form factor \(f(q)\)

Return type:

ndarray

reborn.target.atoms.hubbel_henke_scattering_factors(q_mags, atomic_number, photon_energy=None, beam=None)[source]

Get the q-dependent atomic form factors for a single atomic number and single photon energy, using the Hubbel atomic form factors (Hubbel et al. (1975)) and the Henke dispersion corrections (Henke et al. (1993)).

Parameters:
  • q_mags (ndarray) – q vector magnitudes.

  • atomic_number (int) – Atomic number.

  • photon_energy (float) – Photon energy.

  • beam (Beam) – Supply a beam instead of photon energy.

Returns:

Atomic form factor \(f(q)\) with dispersion corrections

Return type:

ndarray

reborn.target.atoms.cmann_henke_scattering_factors(q_mags, atomic_number, photon_energy=None, beam=None)[source]

Get the q-dependent atomic form factors for a single atomic number and single photon energy, using the Cromer-Mann scattering factors (Cromer et al. (1968)) and the Henke dispersion corrections (Henke et al. (1993)).

Parameters:
  • q_mags (ndarray) – q vector magnitudes.

  • atomic_number (int) – Atomic number.

  • photon_energy (float) – Photon energy.

  • beam (Beam) – Supply a beam instead of photon energy.

Returns:

Atomic form factor \(f(q)\) with dispersion corrections

Return type:

ndarray

reborn.target.atoms.cromer_mann_scattering_factors(q_mags, atomic_numbers)[source]

Get Cromer-Mann scattering factors (Cromer et al. (1968)) for a range of atomic numbers and q magnitudes. The Cromer-Mann formula is

(1)\[\begin{split}f(q) = c + \sum_{i=1}^4 a_i \exp[-b_i q^2 / 16\pi^2] \\\end{split}\]

Note: The Cromer-Mann tables contain the \(a_i\) and \(b_i\) coefficients for ions, but this function only works with neutral atoms for simplicity. You may use get_cmann_form_factors if you want to work with ions.

Parameters:
  • atomic_numbers (ndarray) – Array of atomic numbers. Will be converted to 1D ndarray of int type.

  • q_mags (ndarray) – Array of q magnitudes. Will be converted to 1D ndarray or floats.

Returns:

Atomic scattering factor array with shape (len(atomic_numbers), len(q_mags))

Return type:

ndarray

reborn.target.atoms.get_cmann_form_factors(cman, q)[source]

Generates the form factors \(f_0(q)\) from the model of Cromer et al. (1968). If you want to work with neutral atoms, use cromer_mann_scattering_factors instead.

Parameters:
  • cman (dict) – Dictionary containing {atomic_number : array of cman parameters}

  • q (ndarray) – Either an array of vectors \(\vec{q}\) with shape (N,3), or an array of \(q\) magnitudes.

Returns:

Ask Derek Mendez.

reborn.target.atoms.get_cromermann_parameters(atomic_numbers=None)[source]

Get Cromer-Mann parameters for each atom type.

Modified from tjlane/thor.git by Derek Mendez.

Parameters:

atomic_numbers (ndarray, int) – A numpy array of the atomic numbers of each atom in the system.

Returns:

The Cromer-Mann parameters for the system. The dictionary key corresponds to the atomic number, and the

parameters are listed in order as a_1, a_2, a_3, a_4, b_1, b_2, b_3, b_4, c

Return type:

dict

reborn.target.atoms.get_cromer_mann_densities(Z, r)[source]

Given the Cromer-Mann scattering factors \(f(q)\), we calculate \(\rho(r)\) as follows:

(2)\[\begin{split}\rho(r) &= \frac{1}{8\pi} \int f(q) \exp(i \vec{q}\cdot\vec{r}) d^3 q \\ &= \frac{1}{8\pi} \int \left[ c + \sum_{i=1}^4 a_i \exp(-b_i q^2 / 16\pi^2) \right] \exp(i \vec{ q}\cdot\vec{r}) d^3 q\end{split}\]
Parameters:
  • Z

  • r

Returns:

reborn.target.crystal module

Basic utilities for dealing with crystalline objects.

reborn.target.crystal.get_pdb_file(pdb_id, save_path='/tmp/reborn', silent=False, use_cache=True)[source]

Download a PDB file from the PDB web server and return the path to the downloaded file. There is a data directory included with reborn that includes a few PDB files for testing purposes - if the requested PDB file exists there, the path to the file included with reborn will be returned.

Note: After you download a file for the first time, the downloaded file path will be returned on the next call to this function to avoid downloading the same file multiple times.

Parameters:
  • pdb_id (str) – For example: “101M” or “101M.pdb”. There are also some special strings which so far include “lysozyme” (2LYZ) and “PSI” (1jb0)

  • save_path (string) – Path to the downloaded file. The default is the current working directory, which depends on where you run your code.

  • silent (bool) – Set to true if you do not want any text printed to stdout.

  • use_cache (bool) – Try to find file in reborn’s cache / standard data. Default: True.

Returns:

Path to PDB file

Return type:

string

class reborn.target.crystal.UnitCell(a, b, c, alpha, beta, gamma)[source]

Bases: object

Simple class for primitive unit cell information. Provides the convenience methods r2x(), x2r(), q2h() and h2q() for transforming between fractional and orthogonal coordinates in real space and reciprocal space.

Note

We do not handle non-primitive unit cells in reborn. We always assume that the unit cell contains one and only one crystal lattice point per cell.

Always initialize with the lattice parameters. Units are SI and radians.

Parameters:
  • a – Lattice constant

  • b – Lattice constant

  • c – Lattice constant

  • alpha – Lattice angle

  • beta – Lattice angle

  • gamma – Lattice angle

a = None

Lattice constant (float)

b = None

Lattice constant (float)

c = None

Lattice constant (float)

alpha = None

Lattice angle in radians (float)

beta = None

Lattice angle in radians (float)

gamma = None

Lattice angle in radians (float)

property o_mat

Orthogonalization matrix (3x3 array). Does the transform r = O.x on fractional coordinates x.

property o_mat_inv

Inverse orthogonalization matrix (3x3 array)

property a_mat

Orthogonalization matrix transpose (3x3 array). Does the transform q = A.h, with Miller indices h.

property a_mat_inv

Orthogonalization matrix transpose (3x3 array). Does the transform q = A.h, with Miller indices h.

property volume

Unit cell volume (float)

r2x(r_vecs)[source]

Transform orthogonal coordinates to fractional coordinates.

x2r(x_vecs)[source]

Transform fractional coordinates to orthogonal coordinates.

q2h(q_vecs)[source]

Transform reciprocal coordinates to fractional Miller coordinates.

h2q(h_vecs)[source]

Transform fractional Miller coordinates to reciprocal coordinates.

property a_vec

Crystal basis vector a.

property b_vec

Crystal basis vector b.

property c_vec

Crystal basis vector c.

property as_vec

Reciprocal basis vector a*.

property bs_vec

Reciprocal basis vector b*.

property cs_vec

Reciprocal basis vector c*.

class reborn.target.crystal.SpaceGroup(spacegroup_symbol=None, sym_rotations=None, sym_translations=None, hermann_mauguin_symbol=None, hall_number=None, itoc_number=None)[source]

Bases: object

Container for crystallographic spacegroup information. Most importantly, transformation matices and vectors. These transformations are purely in the fractional coordinate basis. Note that this class has no awareness of the meanings of spacegroup symbols – I have not yet found a good way to programmatically go from a spacegroup symbol string to a set of symmetry operators.

spacegroup_symbol = None

Spacegroup symbol (free format… has no effect)

sym_rotations = None
sym_translations = None
hermann_mauguin_symbol = None
property n_molecules

The number of symmetry operations.

property n_operations

The number of symmetry operations.

apply_symmetry_operation(op_num, x_vecs)[source]

Apply a symmetry operation to an asymmetric unit.

Parameters:
  • op_num (int) – The number of the operation to be applied.

  • x_vecs (Nx3 array) – The atomic coordinates in the crystal basis.

Returns: Nx3 array.

apply_inverse_symmetry_operation(op_num, x_vecs)[source]

Apply an inverse symmetry operation to an asymmetric unit.

Parameters:
  • op_num (int) – The number of the operation to be applied.

  • x_vecs (Nx3 array) – The atomic coordinates in the crystal basis.

Returns: Nx3 array.

class reborn.target.crystal.CrystalStructure(pdb_file_path, no_warnings=False, expand_ncs_coordinates=False, create_bio_assembly=0, tight_packing=False, unitcell=None, spacegroup=None, tempdir=None, ignore_waters=False)[source]

Bases: object

A container class for stuff needed when dealing with crystal structures. Combines information regarding the asymmetric unit, unit cell, and spacegroup.

This class is initialized with a PDB file.

FIXME: Allow initialization without a PDB file.

Parameters:
  • pdb_file_path (string) – Path to a pdb file (or just PDB ID if you want it auto-downloaded)

  • no_warnings (bool) – Suppress warnings concerning ambiguities in symmetry operations

  • expand_ncs_coordinates (bool) – Duplicate asymmetric unit to include non-crystallographic symmetry (NCS) partners.

  • create_bio_assembly (int) – Duplicate asymmetric unit to include biological assembly model. The integer specifies the desired biological assembly to create. See the BIOMOLECULE label of the REMARK section of the PDB (accurate as of Jun 2022)

  • tight_packing (bool) – Choose whether or not to enable a physical arrangement of the asymmetric unit symmetry partners such that each molecular center of mass fits within the unit cell. Useful for finite crystals.

  • unitcell (Unit Cell class) – Specify a unit cell manually if needed.

  • spacegroup (Space Group class) – Specify a space group manually if needed.

  • tempdir (str) – Where to save downloaded pdb file.

  • ignore_waters (bool) – Ignore water molecules, default is False.

mosaicity_fwhm = 0
crystal_size = 1e-06
crystal_size_fwhm = 0.0
mosaic_domain_size = 5e-07
mosaic_domain_size_fwhm = 0.0
cryst1 = ''
pdb_dict = None
unitcell = None

UnitCell class instance

spacegroup = None

Spacegroup class instance

fractional_coordinates = None
molecule = None

Molecule class instance containing the asymmetric unit

property x_vecs

Fractional coordinates of atoms.

Returns: Nx3 numpy array

property x
get_symmetry_expanded_coordinates()[source]

Returns: All atomic coordinates including spacegroup symmetry partners.

set_tight_packing()[source]

Re-define the atomic coordinates and symmetry transforms so that the center-of-mass of each molecule lies within the unit cell. Two steps: first shift the coordinates so that the COM of the asymmetric unit lies in the unit cell. Next, add integral shifts to each of the symmetry operations as needed to put the COM of all symmetry partners into the unit cell.

Returns: None

class reborn.target.crystal.FiniteLattice(max_size=None, unitcell=None)[source]

Bases: object

A utility for creating finite crystal lattices. Uses an occupancy model in which a maximum-size array of occupancies is created and set to 1 (occupied), and subsequently facets may be added by specifying the plane at which the cut is made. A cut sets the occupancies beyond the plane to 0 (unoccupied). Lattice vector positions can then be generated in the crystal or cartesian basis. Gaussian disorder may be added. Special shapes are supported, including hexagonal prisms, parallelepipeds, and spheres.

Parameters:
  • max_size – Integer N that sets the size of the lattice to N x N x N.

  • unitcell – A crystal.UnitCell type that is needed to generate

property occupied_x_coordinates

The occupied coordinate vectors in the crystal basis. An (M, 3) numpy array.

property occupied_r_coordinates

The occupied coordinate vectors in the cartesian (laboratory) basis. An (M, 3) numpy array.

add_facet(plane=None, length=None, shift=0)[source]

Creates a crystal facet by zeroing the lattice occupancies for which the following condition is met:

(3)\[(\mathbf{x}+\mathbf{s})\cdot\mathbf{p} > L\]

where \(\mathbf{x}\) are the nominal lattice positions (with one lattice point sitting on the origin), \(\mathbf{s}\) is an optional shift (which might correspond to the center of mass of one molecule symmetry partner in a spacegroup), \(\mathbf{p}\) defines the facet orientation (a vector that points from the origin to the normal of the facet surface), and \(L\) is the length from the origin to the facet surface.

This operation will set the occupancies of the rejected lattice sites to zero. When you subsequently access the occupied_x_coordinates property, the zero-valued occupancies will not be returned.

Parameters:
  • plane (numpy array) – The vector \(\vec{p}\) defined above.

  • length (numpy array) – The length \(L\) defined above.

  • shift (numpy array) – The vector \(\vec{s}\) defined above.

reset_occupancies()[source]

Set all occupancies to 1

make_hexagonal_prism(width=3, length=10, shift=0)[source]

Specialized function to form a hexagonal prism by adding eight facets. A crude illustration is shown below. This method assumes a “standard” hexagonal cell in which alpha=90, beta=90, gamma=120. The length and width parameters span from facet to facet. The length facets are in the planes [0,0,1] and [0,0,-1]. The three widths specify the facet pairs ([1,0,0], [-1,0,0]), ([0,1,0], [0,-1,0]), ([1,-1,0], [-1,1,0]). Note that, by default, there is always at minimum one lattice point that lies on the origin; if that is not desired then you may use the shift parameter as discussed in the add_facet method.

Parameters:
  • width (float or array) – Three widths to specify the prism shape/size, as explained above.

  • length (float) – Length of the prism, as illustrated above.

  • shift (numpy array) – An optional shift of the entire crystal (defaults to 0).

make_parallelepiped(shape=(5, 5, 5), shift=0)[source]

Cuts out a Parallelepiped shaped crystal.

Parameters:
  • shape (3-tuple) – Number of lattice points to include in the parallelepiped in each direction.

  • shift (numpy array) – An optional shift of the entire crystal (defaults to 0).

make_sphere(radius)[source]

Create a spherical crystal, using the radius in cartesian space.

Parameters:

radius (float) – The cartesian-space radius.

set_gaussian_disorder(sigmas=(0.0, 0.0, 0.0))[source]

Add Gaussian-distributed random offsets to all lattice points

Parameters:

sigmas (numpy array) – The three standard deviations along crystal basis vectors (so these should probably be less than 1)

class reborn.target.crystal.CrystalDensityMap(cryst, resolution, oversampling)[source]

Bases: object

A helper class for working with 3D crystal density maps. Most importantly, it helps with spacegroup symmetry transformations. Once initialized with a crystal spacegroup and lattice, along with desired resolution (\(1/d\)) and oversampling ratio (equal to 1 for normal crystallographic Bragg sampling), this tool chooses the shape of the map and creates lookup tables for the symmetry transformations. The shape is chosen according to the formula \((N_1, N_2, N_3) = \text{roundup}(a_1/d, a_2/d, a_3/d)\) where the function \(\texttt{roundup}\) takes the next-largest integer multiple of \(P\). The value of \(P\) comes from the point group and is equal to the largest number of rotations about any axis. For example, the spacegroup \(P6_3\) has a point group of 6, and thus \(P=6\).

It is generally assumed that you are working in the crystal basis (fractional coordinates). When working in the crystal basis, the symmetry transforms avoid interpolation artifacts. The shape of the map is chosen specifically to avoid interpolation artifacts.

This class does not maintain the data array as the name might suggest. It provides methods needed to work on the data arrays.

On initialization, you provide a CrystalStructure class instance, along with your desired resolution and oversampling.

Parameters:
  • cryst (CrystalStructure instance) – A crystal structure that contains the spacegroup and lattice information.

  • resolution (float) – The desired resolution of the map. The resolution \(d\) is defined such that the number of samples along the \(n`__th edge is greater than or equal to :math:`a_n/d\) for lattice constant \(a_n\). The actual value will be rounded up toward the nearest multiple of 1, 2, 3, 4, or 6 (whichever is appropriate for the spacegroup symmetry).

  • oversampling (int) – An oversampling of 2 gives a real-space map that is twice as large as the unit cell. In Fourier space, there will be one sample between Bragg samples. And so on for 3,4,…

sym_luts = None

Symmetry lookup tables – indices that map AU to sym. partner

cryst = None

CrystalStructure class used to initiate the map

oversampling = None

Oversampling ratio

dx = None

Length increments for fractional coordinates

cshape = None

Number of samples along edges of unit cell within density map

shape = None

Number of samples along edge of full density map (includes oversampling)

size = None

Total number of elements in density map (np.prod(self.shape)

strides = None

The stride vector (mostly for internal use)

property n_vecs

Get an Nx3 array of vectors corresponding to the indices of the map voxels. The array looks like this:

[[0, 0, 0], [0, 0, 1], [0, 0, 2], … , [Nx-1,Ny-1,Nz-1]]

Note that it is the third index, which we might associate with “z”, that increments most rapidly.

Returns: ndarray

property x_vecs

Get an Nx3 array that contains the fractional coordinates. For example, if there were four samples per unit cell, the array will look like this:

[[0, 0, 0], [0, 0, 0.25], [0, 0, 0.5], [0, 0, 0.75], [0, 0.25, 0], … , [ 0.75, 0.75, 0.75]]

Returns: ndarray

property x_limits

Depreciated. Do not use this.

property x_min

Vector pointing to the origin corner of the grid. All zeros.

property x_max

Vector pointing to the far corner of the grid (furthest from the origin).

property h_vecs

This provides an \(M \times 3\) array of Fourier-space vectors \(\vec{h}\) that correspond to the real-space vectors \(\vec{x}\) of this density map. These \(\vec{h}\) vectors can be understood as the “fractional Miller indices”, and they are defined in accordance with the numpy FFT convention . The numpy foward FFT is defined in these terms as

(4)\[F_h = \sum_{x=0}^{N-1} f_x \exp(-2 \pi i h x / N) \; , \quad h = 0, 1, 2, \ldots , N-1\]

while the inverse FFT is

(5)\[f_x = \frac{1}{N}\sum_{h=0}^{N-1} F_h \exp(2 \pi i h x / N) \; , \quad x = 0, 1, 2, \ldots , N-1\]

Note that the above expression is defined such that both \(h\) and \(x\) have integer step size. Our CrystalDensityMap class handles oversampling, in which case the integer \(h\) becomes \(h/s\) and the array size \(N\) becomes :math`sN`, where \(s\) is the oversampling ratio.

Returns: ndarray

property h_limits

This is depreciated. Do not use it.

Returns:

property h_min

The lower limits of the \(\vec{h}\) vectors, which correspond to the FFT of this density map.

Returns:

numpy array

property h_max

The upper limits of the \(\vec{h}\) vectors, which correspond to the FFT of this density map.

Returns:

numpy array

property q_vecs

For convenience, converts self.h_vecs to q vectors (using self.cryst.unitcell.h2q method).

property q_min

For convenience, converts self.h_min to a q vector (using self.cryst.unitcell.h2q method). Note that sampling in q space will not match with the sampling in h space unless the unit cell is orthorhombic.

property q_max

For convenience, converts self.h_max to a q vector (using self.cryst.unitcell.h2q method).

get_sym_luts()[source]

This provides a list of “symmetry transform lookup tables” (symmetry LUTs). These are the flattened array indices that map voxels from the asymmetric unit (AU) map to the kth symmetry partner. This kind of transformation, from AU to symmetry partner k, is performed as follows:

# data is a 3D array of densities
data_trans = np.empty_like(data)
lut = crystal_density_map.get_sym_luts()[k]
data_trans.flat[lut] = data.flat[:]

For convenience the above operation may be performed by the method au_to_k, while the inverse operation may be performed by the method k_to_au. The method symmetry_transform can be used to transform from one symmetry partner to another.

Note that the LUTs are kept in memory for future use - beware of the memory requirement.

Returns:

The symmetry lookup tables (LUTs)

Return type:

list of numpy arrays

au_to_k(k, data)[source]

Transform a map of the asymmetric unit (AU) to the kth symmetry partner. Note that the generation of the first symmetry partner (k=0, where k = 0, 1, …, N-1) might employ a non-identity rotation matrix and/or a non-zero translation vector – typically this is not the case but it can happen for example if the symmetry operations are chosen such that all molecules are packed within the unit cell.

Parameters:
  • k (int) – The index of the symmetry partner (starting with k=0)

  • data (ndarray) – The input data array.

Returns:

Transformed array

Return type:

ndarray

k_to_au(k, data)[source]

Transform a map of the kth symmetry partner to the asymmetric unit. This reverses the action of the au_to_k method.

Parameters:
  • k (int) – The index of the symmetry partner (starting with k=0)

  • data (3D ndarray) – The input data array.

Returns:

Transformed array

Return type:

3D ndarray

symmetry_transform(i, j, data)[source]

Apply crystallographic symmetry transformation to a density map (3D ndarray). This applies the mapping from symmetry element i to symmetry element j, where i=0,1,…,N-1 for a spacegroup with N symmetry operations.

Parameters:
  • i (int) – The “from” index; symmetry transforms are performed from this index to the j index

  • j (int) – The “to” index; symmetry transforms are performed from the i index to this index

  • data (ndarray) – The 3D block of data to apply the symmetry transform on

Returns: ndarray with transformed densities

place_atoms_in_map(atom_x_vecs, atom_fs, mode='trilinear', fixed_atom_sigma=5e-11)[source]

This will take a list of atom position vectors and densities and place them in a 3D map. The position vectors should be in the crystal basis, and the densities must be real.

FIXME: This function will soon be replaced with one that utilizes atomic form factors.

Parameters:
  • atom_x_vecs (numpy array) – An nx3 array of position vectors

  • atom_fs (numpy array) – An n-length array of densities (must be real)

  • mode (str) – Either ‘nearest’, ‘trilinear’, or ‘gaussian’ (default: ‘trilinear’)

  • fixed_atom_sigma (float) – Standard deviation of the Gaussian atoms

Returns:

The sum of densities that were provided as input.

Return type:

numpy array

property voxel_volume

Get the volume of a voxel in orthogonal coordinates (i.e. Cartesian).

Returns:

Voxel volume

Return type:

float

reborn.target.crystal.build_atomic_scattering_density_map(x_vecs, f, sigma, x_min, x_max, shape, orth_mat, n_sigma)[source]
reborn.target.crystal.place_atoms_in_map(x_vecs, atom_fs, sigma, s, orth_mat, map_x_vecs, f_map, f_map_tmp)[source]

Needs documentation…

reborn.target.crystal.pdb_to_dict(pdb_file_path, ignore_waters=False, verbose=False)[source]

Return a dictionary with a subset of PDB information. If there are multiple atomic models, only the first will be extracted. Units are the standard PDB units: angstrom and degrees.

Parameters:
  • pdb_file_path – Path to pdb file

  • ignore_waters – Ignore water molecules if True. (Default: False)

Returns:

Rick: Think about better formatting of returns. A dictionary with the following keys: - ‘scale_matrix’ (from SCALE) - ‘scale_translation’ (from SCALE) - ‘atomic_coordinates’ (from ATOM or HETATM) - ‘atomic_symbols’ (from ATOM or HETATM) - ‘unit_cell’ (from CRYST1) - ‘spacegroup_symbol’ (from CRYST1) - ‘spacegroup_rotations’ (from SMTRY) - ‘spacegroup_translations’ (from SMTRY) - ‘ncs_rotations’ (from MTRIX) - ‘ncs_translations’ (from MTRIX)

Return type:

FIXME

class reborn.target.crystal.FiniteCrystal(cryst, max_size=20)[source]

Bases: object

Utility that allows for the shaping of finite crystal lattices, with consideration of the crystal structure (molecular structure, spacegroup, etc.). This is useful for simulating complete crystals with strange edges or other defects that depart from idealized crystals.

Parameters:
  • cryst (CrystalStructure) – A crystal structure object. The center-of-mass of asymmetric unit and spacegroup provided by this object will affect the centering of the lattices.

  • max_size (3-element array) – Same as in the FiniteLattice class.

cryst = None
au_x_coms = None
lattices = None
add_facet(plane=None, length=None)[source]

See equivalent method in FiniteLattice. In this case, the facet is added to all of the finite lattices (one for each symmetry partner).

reset_occupancies()[source]

See equivalent method in FiniteLattice. In this case, all occupancies are reset (i.e. for all symemtry partner lattices).

make_hexagonal_prism(width=3, length=10)[source]

See equivalent method in FiniteLattice. In this case, the facets are added to all of the finite lattices (one for each symmetry partner).

make_parallelepiped(shape=(5, 5, 5), shift=0)[source]

See equivalent method in FiniteLattice. In this case, the facets are added to all of the finite lattices (one for each symmetry partner).

set_gaussian_disorder(sigmas=(0, 0, 0))[source]

See equivalent method in FiniteLattice .

reborn.target.crystal.get_hm_symbols()[source]

Returns a list of the 530 Hermann Mauguin symbols.

reborn.target.crystal.hermann_mauguin_symbol_from_hall_number(hall_number)[source]
reborn.target.crystal.itoc_number_from_hall_number(hall_number)[source]
reborn.target.crystal.hall_number_from_hermann_mauguin_symbol(hermann_mauguin_symbol)[source]
reborn.target.crystal.itoc_number_from_hermann_mauguin_symbol(hermann_mauguin_symbol)[source]
reborn.target.crystal.hall_number_from_itoc_number(itoc_number)[source]
reborn.target.crystal.hermann_mauguin_symbol_from_itoc_number(itoc_number)[source]
reborn.target.crystal.spacegroup_ops_from_hall_number(hall_number)[source]
reborn.target.crystal.get_symmetry_operators_from_hall_number(hall_number)[source]
reborn.target.crystal.get_symmetry_operators_from_space_group(hm_symbol)[source]

For a given Hermann-Mauguin spacegroup, provide the symmetry operators in the form of translation and rotation operators. These operators are in the crystal basis. As of now, this is a rather unreliable function. My brain hurts badly after attempting to make sense of the ways in which space groups are specified… if you are having trouble finding your spacegroup check that it is contained in the output of print(crystal.get_hm_symbols()).

Input:

hm_symbol (str): A string indicating the spacegroup in Hermann–Mauguin notation (as found in a PDB file)

Output: Two lists: Rs and Ts. These correspond to lists of rotation matrices (3x3 numpy

arrays) and translation vectors (1x3 numpy arrays)

reborn.target.density module

class reborn.target.density.DensityMap(shape=None, corner_min=None, corner_max=None, deltas=None)[source]

Bases: object

Parameters:
  • shape – shape of the 3d map

  • corner_min – a 3-vector with the coordinate of the center of the corner pixel

  • corner_max – opposit corner to corner_min

  • deltas – step size along each of the 3 dimensions

spacings = None
corner_min = None
shape = None
corner_max = None
n_voxels = None
dx = None
strides = None
property n_vecs

Get an Nx3 array of vectors corresponding to the indices of the map voxels. The array looks like this:

[[0,0,0],[0,0,1],[0,0,2], … ,[N-1,N-1,N-1]]

Returns: ndarray

reborn.target.density.build_atomic_scattering_density_map(x_vecs, f, sigma, x_min, x_max, shape, orth_mat, max_radius=10000000000.0)[source]

Construct an atomic scattering density by summing Gaussians. Sampling is assumed to be a rectangular grid, but axes need not be orthogonal (the orthogonalization matrix may be provided). Normalization is preformed such that the sum over the whole map should equal the sum over the input scattering factors f.

In order to increase speed, you may specify the maximum atom size, in which case only grid points within that distance are sampled.

Parameters:
  • x_vecs (Mx3 numpy array) – Atom position vectors

  • f (numpy array) – Atom structure factors

  • sigma (float) – Standard deviation of Gaussian (all atoms are the same)

  • x_min (numpy array) – Min position of corner voxel (center of voxel)

  • x_max (numpy array) – Max position of corner voxel (center of voxel)

  • shape (numpy array) – Shape of the 3D array

  • orth_mat (3x3 numpy array) – Matrix that acts on distance vectors before calculating distance scalar

  • max_atom_size – Maximum atom size (saves time by not computing tails of Gaussians)

Returns:

The density map

Return type:

numpy array

reborn.target.density.trilinear_interpolation(*args, **kwargs)[source]

Depreciated duplicate of reborn.misc.interpolate.trilinear_interpolation

reborn.target.density.trilinear_insertion(*args, **kwargs)[source]

Depreciated duplicate of reborn.misc.interpolate.trilinear_insertion

reborn.target.density.trilinear_insertions(*args, **kwargs)[source]

Depreciated duplicate of reborn.misc.interpolate.trilinear_insertions

reborn.target.density.trilinear_insertion_factor(*args, **kwargs)[source]

Depreciated duplicate of reborn.misc.interpolate.trilinear_insertion_factor

reborn.target.molecule module

Utilities for manipulating molecules.

class reborn.target.molecule.Molecule(coordinates=None, atomic_numbers=None, atomic_symbols=None)[source]

Bases: object

Simple container for organizing atomic numbers, coordinates, etc. There is presently no information about bonds, as the name of the class might suggest…

Parameters:
  • coordinates (numpy array) – An array of shape Nx3 vectors in cartesiaon coordinates

  • atomic_numbers (numpy array) – Array of integer atomic numbers

  • atomic_symbols (numpy array) – Array of atomic symbols (as an alternative to providing atomic numbers)

atomic_numbers = None
n_atoms = None
occupancies = None
property coordinates
get_scattering_factors(photon_energy=None, beam=None)[source]

Get atomic scattering factors. You need to specify the photon energy or pass a Beam class instance in.

This wraps the function reborn.target.atoms.get_scattering_factors for more details; see the docs for more details.

Parameters:
  • photon_energy – In SI units as always

  • beam – Optionally, provide a reborn.beam.Beam instance instead of photon_energy

Returns:

Numpy array of complex scattering factors

property max_atomic_pair_distance

Return the maximum distance between two atoms in the molecule. This was intended to be useful for determining an appropriate angular or q-space sampling frequency.

get_size()[source]

Quickly get molecule size (max edge of enclosing rectangle).

get_molecular_weight()[source]

Returns the molecular weight in SI units (kg).

get_centered_coordinates()[source]

Get the coordinates with center of mass set to the origin.

Returns:

An Nx3 array of coordinates.

Return type:

ndarray

get_atomic_weights()[source]

Returns the atomic weights in SI units (kg).

view()[source]
remove_hydrogen_atoms()[source]
reborn.target.molecule.get_molecule(molecule_type=None)[source]

Convenience function to create common molecules. Currently, the known types are: N2, O2, H2O, or any of the elements in the periodic table.

Parameters:

molecule_type (str, int, Molecule) – Molecule type

Returns:

Molecule

reborn.target.placer module

Utilities for positioning objects (x-ray targets) in space.

class reborn.target.placer.Place(box_edge, min_dist, max_try=10000, *args, **kwargs)[source]

Bases: cKDTree

Utility for placing points into a box while ensuring that no two points get closer than some minimum distance. This is based on the scipy.spatial.cKDTree class.

Note

  • Since an iterative process is employed, make sure that you use reasonable values for the minimum distance and box size. It is of course possible to choose parameters for which it is impossible to place the spheres.

  • Since this is a sub-class of scipy.spatial.cKDTree, you may pass arguments that are relevant to that class.

Contributed by Derek Mendez.

Parameters:
  • box_edge (float) – Side length of the box to place spheres into.

  • min_dist (float) – Minimum distance between two points in the box.

  • max_try (int) – Number of times to try placing a new point such that is does not overlap.

insert()[source]

Adds one new point to the box.

reborn.target.placer.place_spheres(volume_fraction, radius=1.0, box_edge=None, n_spheres=1000, tol=0.01)[source]

Contributed by Derek Mendez.

No documentation.

Parameters:
  • volume_fraction (float) – Fraction of sample volume occupied by spheres.

  • radius (float) – Undocumented.

  • box_edge (float) – Undocumented.

  • n_spheres (int) – How many spheres in the sample volume.

  • tol (float) – Minimum distance the unit spheres can be to one another.

Returns: Ask Derek.

reborn.target.placer.particles_in_a_sphere(sphere_diameter, n_particles, particle_diameter, max_attempts=1000000.0)[source]

Place particles randomly in a spherical volume. Assume particles are spheres and they cannot touch any other particle. Also assumes that surface of spherical particles cannot extend beyond the surface of the containing sphere (thus, the maximum distance of a particle from the origin is sphere_diameter/2 - particle_diameter/2).

Parameters:
  • sphere_diameter (float) – Diameter of the bounding sphere, within which particle positions must lie (but without the particle surface extending beyond this bounding sphere).

  • n_particles (int) – Number of particles to fit in the bounding sphere.

  • particle_diameter (float) – Diameter of the particles that must fit in the sphere.

  • max_attempts (int) – Optional. How many times to try to place a sphere before giving up. Default: 1e6.

Returns:

The array of vectors with particle positions.

Return type:

ndarray

reborn.target.placer.particles_in_a_cylinder(cylinder_diameter, cylinder_length, n_particles, particle_diameter, max_attempts=1000000.0)[source]

Place particles randomly in a cylindrical volume. Assume particles are spheres and they cannot touch any other particle. Also assumes that surface of spherical particles cannot extend beyond the surface of the containing cylinder (thus, the maximum distance of a particle from the axis is cylinder_diameter/2 - particle_diameter/2).

Parameters:
  • cylinder_diameter (float) – Diameter of the bounding sphere, within which particle positions must lie (but without the particle surface extending beyond this bounding sphere).

  • n_particles (int) – Number of particles to fit in the bounding sphere.

  • particle_diameter (float) – Diameter of the particles that must fit in the sphere.

  • max_attempts (int) – Optional. How many times to try to place a sphere before giving up. Default: 1e6.

Returns:

The array of vectors with particle positions.

Return type:

ndarray

Module contents

The target package contains classes, functions, etc. that are related to the objects we shoot with x-rays. Crystals, molecules, atoms, electron density maps, and so on.