Saxs#

class maicos.Saxs(atomgroup: AtomGroup, unwrap: bool = False, refgroup: AtomGroup | None = None, jitter: float = 0.0, concfreq: int = 0, bin_spectrum: bool = True, qmin: float = 0, qmax: float = 6, dq: float = 0.1, thetamin: float = 0, thetamax: float = 180, output: str = 'sq.dat')[source]#

Bases: AnalysisBase

Small angle X-Ray scattering intensities (SAXS).

This module computes the structure factor \(S(q)\), the scattering intensity \(I(q)\) and their corresponding scattering vectors \(q\). For a system containing only one element the structure factor and the scattering intensity are connected via the form factor \(f(q)\)

\[I(q) = [f(q)]^2 S(q)\]

For more details on the theory behind this module see Small-angle X-ray scattering.

By default the scattering vectors \(\boldsymbol{q}\) are binned according to their length \(q\) using a bin width given by dq. Setting the option bin_spectrum=False, also the raw scattering vectors and their corresponding Miller indices can be saved. Saving the scattering vectors and Miller indices is only possible when the box vectors are constant in the whole trajectory (NVT) since for changing cells the same Miller indices correspond to different scattering vectors.

Analyzed scattering vectors \(q\) can be restricted by a minimal and maximal angle with the z-axis. For 0 and 180, all possible vectors are taken into account. To obtain the scattering intensities, the structure factor is normalized by an element-specific form factor based on Cromer-Mann parameters Prince[1].

For the correlation time estimation the module will use the value of the scattering intensity with the largest possible \(q\) value.

For an example on the usage refer to How-to: SAXS.

Parameters:
  • atomgroup (MDAnalysis.core.groups.AtomGroup) – A AtomGroup for which the calculations are performed.

  • unwrap (bool) –

    When True, molecules that are broken due to the periodic boundary conditions are made whole.

    If the input contains molecules that are already whole, speed up the calculation by disabling unwrap. To do so, use the flag -no-unwrap when using MAICoS from the command line, or use unwrap=False when using MAICoS from the Python interpreter.

    Note: Molecules containing virtual sites (e.g. TIP4P water models) are not currently supported in MDAnalysis. In this case, you need to provide unwrapped trajectory files directly, and disable unwrap. Trajectories can be unwrapped, for example, using the trjconv command of GROMACS.

  • refgroup (MDAnalysis.core.groups.AtomGroup) – Reference AtomGroup used for the calculation. If refgroup is provided, the calculation is performed relative to the center of mass of the AtomGroup. If refgroup is None the calculations are performed with respect to the center of the (changing) box.

  • jitter (float) –

    Magnitude of the random noise to add to the atomic positions.

    A jitter can be used to stabilize the aliasing effects sometimes appearing when histogramming data. The jitter value should be about the precision of the trajectory. In that case, using jitter will not alter the results of the histogram. If jitter = 0.0 (default), the original atomic positions are kept unchanged.

    You can estimate the precision of the positions in your trajectory with maicos.lib.util.trajectory_precision(). Note that if the precision is not the same for all frames, the smallest precision should be used.

  • concfreq (int) – When concfreq (for conclude frequency) is larger than 0, the conclude function is called and the output files are written every concfreq frames.

  • bin_spectrum (bool) – Bin the spectrum. If False Miller indices of q-vector are returned. Only works for NVT simulations.

  • qmin (float) – Starting q (1/Å)

  • qmax (float) – Ending q (1/Å)

  • dq (float) – bin_width (1/Å)

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

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

  • output (str) – Output filename.

results.scattering_vectors#

Length of the binned scattering vectors.

Type:

numpy.ndarray

results.miller_indices#

Miller indices of q-vector (only available if bin_spectrum==False).

Type:

numpy.ndarray

results.struture_factors#

structure factors \(S(q)\)

Type:

numpy.ndarray

results.scattering_intensities#

scattering intensities \(I(q)\)

Type:

numpy.ndarray

References

save() None[source]#

Save results of analysis to file specified by output.