StructureFactor¶
- class mdhelper.analysis.structure.StructureFactor(groups: AtomGroup | tuple[AtomGroup], groupings: str | tuple[str] = 'atoms', *, mode: str = None, form: str = 'exp', dimensions: ndarray[float] | Quantity | Quantity = None, n_points: int = 32, n_surfaces: int = None, n_surface_points: int = 8, q_max: float | Quantity | Quantity = None, wavevectors: ndarray[float] = None, sort: bool = True, unique: bool = True, parallel: bool = False, verbose: bool = True, **kwargs)[source]¶
Bases:
NumbaAnalysisBase
Serial and parallel implementations to calculate the static structure factor \(S(q)\) or partial structure factor \(S_{\alpha\beta}(q)\) for species \(\alpha\) and \(\beta\).
The static structure factor is a measure of how a material scatters incident radiation, and can be computed directly from a molecular dynamics simulation trajectory using
\[\begin{split}S(q)&=\frac{1}{N}\left\langle\sum_{j=1}^N\sum_{k=1}^N \exp{[-i\mathbf{q}\cdot(\mathbf{r}_j-\mathbf{r}_k)]} \right\rangle\\&=\frac{1}{N}\left\langle\left( \sum_{j=1}^N\sin{(\mathbf{q}\cdot\mathbf{r}_j)}\right)^2 +\left(\sum_{j=1}^N\cos{(\mathbf{q}\cdot\mathbf{r}_j)} \right)^2\right\rangle\end{split}\]where \(N\) is the total number of particles or centers of mass, \(\mathbf{q}\) is the scattering wavevector, and \(\mathbf{r}_i\) is the position of the \(i\)-th particle.
For multicomponent systems, the equation above can be generalized to get the partial structure factors
\[\begin{split}S_{\alpha\beta}(q)&=\frac{2-\delta_{\alpha\beta}}{N} \left\langle\sum_{j=1}^{N_\alpha}\sum_{k=1}^{N_\beta} \exp{[-i\mathbf{q}\cdot(\mathbf{r}_j-\mathbf{r}_k)]} \right\rangle\\ &=\frac{2-\delta_{\alpha\beta}}{N}\left\langle \sum_{j=1}^{N_\alpha}\cos{(\mathbf{q}\cdot\mathbf{r}_j)} \sum_{k=1}^{N_\beta}\cos{(\mathbf{q}\cdot\mathbf{r}_k)} +\sum_{j=1}^{N_\alpha}\sin{(\mathbf{q}\cdot\mathbf{r}_j)} \sum_{k=1}^{N_\beta}\sin{(\mathbf{q}\cdot\mathbf{r}_k)} \right\rangle\end{split}\]where \(\delta_{ij}\) is the Kronecker delta, \(N_\alpha\) and \(N_\beta\) are the numbers of particles or centers of mass for species \(\alpha\) and \(\beta\).
The partial structure factors \(S_{\alpha\beta}(q)\) and the static structure factor \(S(q)\) are related via
\[S(q)=\sum_{\alpha=1}^{N_\mathrm{g}} \sum_{\beta=\alpha}^{N_\mathrm{g}}S_{\alpha\beta}(q)\]- Parameters:
- groupsMDAnalysis.AtomGroup or array-like
Group(s) of atoms that share the same grouping type. If
mode=None
, all atoms in the universe must be assigned to a group. Ifmode="pair"
, there must be exactly one or two groups.- groupingsstr or array-like, default:
"atoms"
Determines whether the centers of mass are used in lieu of individual atom positions. If groupings is a str, the same value is used for all groups.
Note
If the desired grouping is not
"atoms"
,the trajectory file should have segments (or chains) containing residues (or molecules) and residues containing atoms, and
residues and segments should be locally unwrapped at the simulation box edges, if not already, using
MDAnalysis.transformations.wrap.unwrap
,MDAnalysis.core.groups.AtomGroup.unwrap()
, orMDAnalysis.lib.mdamath.make_whole()
.
Valid values:
"atoms"
: Atom positions (generally or for coarse-grained simulations)."residues"
: Residues’ centers of mass (for atomistic simulations)."segments"
: Segments’ centers of mass (for atomistic polymer simulations).
- modestr, keyword-only, optional
Evaluation mode.
Valid values:
None
: The static structure factor is computed."pair"
: The partial structure factor is computed between the group(s) in groups."partial"
: The partial structure factors for all unique pairs in groups is computed.
- formstr, keyword-only, default:
"exp"
Expression used to evaluate the structure factors.
Valid values:
"exp"
: Exponential form. Slightly faster due to fewer mathematical operations."trig"
: Trigonometric form. Slightly slower but doesn’t have overflow issues.
- dimensionsarray-like, openmm.unit.Quantity, or pint.Quantity, keyword-only, optional
System dimensions. If not provided, they are retrieved from the topology or trajectory. Only necessary if wavevectors is not specified.
Shape: \((3,)\).
Reference unit: \(\mathrm{Å}\).
- n_pointsint, keyword-only, default:
32
Number of points in the scattering wavevector grid. Additional wavevectors can be introduced via n_surfaces and n_surface_points for more accurate structure factors at small wavenumbers. Alternatively, the desired wavevectors can be specified directly in wavevectors.
- n_surfacesint, keyword-only, optional
Number of spherical surfaces in the first octant that intersect with the grid wavevectors along the three coordinate axes for which to introduce extra wavevectors for more accurate structure factor values. Only available if the system is perfectly cubic.
- n_surface_pointsint, keyword-only, default:
8
Number of extra wavevectors to introduce per spherical surface. Has no effect if n_surfaces is not specified.
- q_maxfloat, openmm.unit.Quantity, or pint.Quantity, keyword-only, optional
Maximum scattering wavevector magnitude.
Reference unit: \(\mathrm{Å}^{-1}\).
- wavevectorsnumpy.ndarray, keyword-only, optional
Scattering wavevectors for which to compute structure factors. Has precedence over n_points, n_surfaces, and n_surface_points if specified.
Reference unit: \(\mathrm{Å}^{-1}\).
- sortbool, keyword-only, default:
True
Determines whether the results are sorted by the wavenumbers.
- uniquebool, keyword-only, default:
True
Determines whether structure factors for the same wavenumber are grouped and averaged.
- parallelbool, keyword-only, default:
False
Determines whether the analysis is performed in parallel.
- verbosebool, keyword-only, default:
True
Determines whether detailed progress is shown.
- **kwargs
Additional keyword arguments to pass to
MDAnalysis.analysis.base.AnalysisBase
.
- Attributes:
- universeMDAnalysis.Universe
MDAnalysis.core.universe.Universe
object containing all information describing the system.- results.unitsdict
Reference units for the results. For example, to get the reference units for
results.wavenumbers
, callresults.units["results.wavenumbers"]
.- results.pairstuple
All unique pairs of indices of the groups of atoms in groups. The ordering coincides with the column indices in results.ssf.
- results.wavenumbersnumpy.ndarray
Wavenumbers \(q\).
Shape: \((N_q,)\).
Reference unit: \(\mathrm{Å}^{-1}\).
- results.ssfnumpy.ndarray
Static structure factor \(S(q)\) or partial structure factors \(S_{\alpha\beta}(q)\).
Shape: \((1,\,N_q)\) or \((C(N_\mathrm{g}+1,\,2),\,N_q)\).
Methods
Computes the partial structure factors given two two-dimensional NumPy arrays, each containing \(\mathbf{q}\cdot\mathbf{r}\), using the trigonometric form.
Performs the calculation.
Saves results to a binary or archive file in NumPy format.
Computes the static structure factors using a two-dimensional NumPy array containing \(\mathbf{q}\cdot\mathbf{r}\) using the trigonometric form.
- static psf_trigonometric_2d_2d(qrs1: ndarray[float], qrs2: ndarray[float]) ndarray[float] [source]¶
Computes the partial structure factors given two two-dimensional NumPy arrays, each containing \(\mathbf{q}\cdot\mathbf{r}\), using the trigonometric form.
\[\frac{NS_{\alpha\beta}(q)}{2-\delta_{\alpha\beta}} =\left\langle \sum_{j=1}^{N_\alpha}\cos{(\mathbf{q}\cdot\mathbf{r}_j)} \sum_{k=1}^{N_\beta}\cos{(\mathbf{q}\cdot\mathbf{r}_k)} +\sum_{j=1}^{N_\alpha}\sin{(\mathbf{q}\cdot\mathbf{r}_j)} \sum_{k=1}^{N_\beta}\sin{(\mathbf{q}\cdot\mathbf{r}_k)} \right\rangle\]- Parameters:
- qrs1np.ndarray
First set of inner products \(\mathbf{q}\cdot\mathbf{r}_j\).
Shape: \((N_q,\,N_r)\).
- qrs2np.ndarray
Second set of inner products \(\mathbf{q}\cdot\mathbf{r}_k\).
Shape: \((N_q,\,N_r)\).
- Returns:
- ssfnp.ndarray
Partial structure factors.
Shape: \((N_q,)\).
- run(start: int = None, stop: int = None, step: int = None, frames: slice | ndarray[int] = None, n_threads: int = None, verbose: bool = None, **kwargs) NumbaAnalysisBase ¶
Performs the calculation.
- Parameters:
- startint, optional
Starting frame for analysis.
- stopint, optional
Ending frame for analysis.
- stepint, optional
Number of frames to skip between each analyzed frame.
- framesslice or array-like, optional
Index or logical array of the desired trajectory frames.
- n_threadsint, keyword-only, optional
Number of threads to use for analysis.
- verbosebool, optional
Determines whether detailed progress is shown.
- **kwargs
Additional keyword arguments to pass to
MDAnalysis.lib.log.ProgressBar
.
- Returns:
- selfNumbaAnalysisBase
Analysis object with results.
- save(file: str | TextIO, archive: bool = True, compress: bool = True, **kwargs) None ¶
Saves results to a binary or archive file in NumPy format.
- Parameters:
- filestr or file
Filename or file-like object where the data will be saved. If file is a str, the
.npy
or.npz
extension will be appended automatically if not already present.- archivebool, default:
True
Determines whether the results are saved to a single archive file. If True, the data is stored in a
.npz
file. Otherwise, the data is saved to multiple.npy
files.- compressbool, default:
True
Determines whether the
.npz
file is compressed. Has no effect whenarchive=False
.- **kwargs
Additional keyword arguments to pass to
numpy.save()
,numpy.savez()
, ornumpy.savez_compressed()
, depending on the values of archive and compress.
- static ssf_trigonometric_2d(qrs: ndarray[float]) ndarray[float] [source]¶
Computes the static structure factors using a two-dimensional NumPy array containing \(\mathbf{q}\cdot\mathbf{r}\) using the trigonometric form.
\[S(q)=\frac{1}{N}\left\langle\left(\sum_{j=1}^N \cos{(\mathbf{q}\cdot\mathbf{r}_j)}\right)^2+\left( \sum_{j=1}^N\sin{(\mathbf{q}\cdot\mathbf{r}_j)} \right)^2\right\rangle\]- Parameters:
- qrsnp.ndarray
Inner products \(\mathbf{q}\cdot\mathbf{r}_j\).
Shape: \((N_q,\,N_r)\).
- Returns:
- ssfnp.ndarray
Static structure factors.
Shape: \((N_q,)\).