Mathematical helper functions

Helper functions for mathematical and physical operations.

scatterkit.lib.math.atomic_form_factor(q: float, element: str) float[source]

Calculate atomic form factor \(f(q)\) for X-ray scattering.

The atomic form factor \(f(q)\) is a measure of the scattering amplitude of a wave by an isolated atom

Attention

The atomic form factor should not be confused with the atomic scattering factor or intensity (often anonymously called form factor). The scattering intensity depends strongly on the distribution of atoms and can be computed using scatterkit.Saxs.

Here, \(f(q)\) is computed in terms of the scattering vector as

\[f(q) = \sum_{i=1}^4 a_i e^{-b_i q^2/(4\pi)^2} + c \,.\]

The coefficients \(a_{1,\dots,4}\), \(b_{1,\dots,4}\) and \(c\) are also known as Cromer-Mann X-ray scattering factors and are documented in Prince[1] and taken from the TU Graz and stored in scatterkit.lib.tables.CM_parameters.

Parameters:
  • q (float) – The magnitude of the scattering vector in reciprocal angstroms (1/Å).

  • element (str) –

    The element for which the atomic form factor is calculated. Known elements are listed in the scatterkit.lib.tables.elements set. United-atom models such as "CH1", "CH2", "CH3", "CH4", "NH1", "NH2", and "NH3" are also supported.

    Note

    element is converted to title case to avoid most common issues with MDAnalysis which uses upper case elements by default. For example "MG" will be converted to "Mg".

Returns:

The calculated atomic form factor for the specified element and q in units of electrons.

Return type:

float

scatterkit.lib.math.compute_rdf_structure_factor(rdf: ndarray, r: ndarray, density: float) tuple[ndarray, ndarray][source]

Computes the structure factor based on the radial distribution function (RDF).

The structure factor \(S(q)\) based on the RDF \(g(r)\) is given by

\[S(q) = 1 + 4 \pi \rho \int_0^\infty \mathrm{d}r r \frac{\sin(qr)}{q} (g(r) - 1)\,\]

where \(q\) is the magnitude of the scattering vector. The calculation is performed via a discrete sine transform as implemented in scipy.fftpack.dst().

For an example take a look at Small-angle X-ray scattering.

Parameters:
  • rdf (numpy.ndarray) – radial distribution function

  • r (numpy.ndarray) – equally spaced distance array on which rdf is defined

  • density (float) – number density of particles

Returns:

  • q (numpy.ndarray) – array of q points

  • struct_factor (numpy.ndarray) – structure factor

scatterkit.lib.math.rdf_structure_factor(rdf: ndarray, r: ndarray, density: float) tuple[ndarray, ndarray][source]

Computes the structure factor based on the radial distribution function (RDF).

The structure factor \(S(q)\) based on an RDF \(g(r)\) is given by

\[S(q) = 1 + 4 \pi \rho \int_0^\infty \mathrm{d}r r \frac{\sin(qr)}{q} (g(r) - 1)\,\]

where \(q\) is the magnitude of the scattering vector. The calculation is performed via a discrete sine transform as implemented in scipy.fftpack.dst().

For an example take a look at Small-angle X-ray scattering.

Parameters:
  • rdf (numpy.ndarray) – radial distribution function

  • r (numpy.ndarray) – equally spaced distance array on which rdf is defined

  • density (float) – number density of particles

Returns:

  • q (numpy.ndarray) – array of q points

  • struct_factor (numpy.ndarray) – structure factor

Raises:

ValueError – If the distance array r is not equally spaced.

scatterkit.lib._cmath.compute_structure_factor(positions, dimensions, qmin, qmax, thetamin, thetamax, weights)

Calculates scattering vectors and corresponding structure factors.

Use via from scatterkit.lib.math import structure_factor

The structure factors are calculated according to

\[S(\boldsymbol{q}) = \left [ \sum\limits_{k=1}^N w_k \cos(\boldsymbol{qr}_k) \right ]^2 + \left [ \sum\limits_{k=1}^N w_k \sin(\boldsymbol{qr}_k) \right ]^2 \,.\]

where \(\boldsymbol{r}_j\) is the positions vector of particle \(k\), \(\boldsymbol{q}\) is scattering vector and the \(w_k\) are optional weights. The possible scattering vectors are determined by the given cell dimensions.

Results are returned as arrays with three dimensions, where the index of each dimensions refers to the Miller indices \(hkl\). Based on the Miller indices and the returned length of the scattering vector the actual scattering vector can be obtained by

\[q_{hkl} = \vert \boldsymbol{q} \vert \frac{2\pi}{L_{hkl}}\]

where \(\vert \boldsymbol{q} \vert\) are the returned lengths of the scattering vector and \(L_{hkl}\) are the components of the simulation cell.

Parameters:
  • positions (numpy.ndarray) – Position array

  • dimensions (numpy.ndarray) – Dimensions of the cell

  • qmin (float) – Starting scattering vector length (1/Å). Possible values are in the range \([0, 180]\).

  • qmax (float) – Ending scattering vector length (1/Å). Possible values are in the range \([0, 180]\).

  • thetamin (float) – Minimal angle (°) between the scattering vectors and the z-axis.

  • thetamax (float) – Maximal angle (°) between the scattering vectors and the z-axis.

  • weights (numpy.ndarray) – Atomic quantity for weighting the structure factor. Provide an array of ones that has the same size as the positions, i.e np.ones(len(positions)), for the standard structure factor.

Returns:

Scattering vectors and their corresponding structure factors.

Return type:

tuple(numpy.ndarray, numpy.ndarray)