DielectricSpectrum

class maicos.DielectricSpectrum(atomgroup: AtomGroup, temperature: float = 300, segs: int = 20, df: float | None = None, bins: int = 200, binafter: float = 20, nobin: bool = False, refgroup: AtomGroup | None = None, unwrap: bool = True, pack: bool = True, jitter: float = 0.0, concfreq: int = 0, output_prefix: str = '')[source]

Bases: AnalysisBase

Linear dielectric spectrum.

This module, given a molecular dynamics trajectory, produces a .txt file containing the complex dielectric function as a function of the (linear, not radial - i.e., \(\nu\) or \(f\), rather than \(\omega\)) frequency, along with the associated standard deviations. The algorithm is based on the Fluctuation Dissipation Relation: \(\chi(f) = -1/(3 V k_B T \varepsilon_0) \mathcal{L}[\theta(t) \langle P(0) dP(t)/dt\rangle]\), where \(\mathcal{L}\) is the Laplace transformation.

Note

The polarization time series and the average system volume are also saved.

Please read and cite [1].

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

  • temperature (float) – Reference temperature (K)

  • segs (int) – Sets the number of segments the trajectory is broken into.

  • df (float) – The desired frequency spacing in THz. This determines the minimum frequency about which there is data. Overrides segs option.

  • bins (int) – Determines the number of bins used for data averaging; (this parameter sets the upper limit). The data are by default binned logarithmically. This helps to reduce noise, particularly in the high-frequency domain, and also prevents plot files from being too large.

  • binafter (int) – The number of low-frequency data points that are left unbinned.

  • nobin (bool) – Prevents the data from being binned altogether. This can result in very large plot files and errors.

  • 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.

  • 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.

  • pack (bool) –

    When True, molecules are put back into the unit cell. This is required because MAICoS only takes into account molecules that are inside the unit cell.

    If the input contains molecules that are already packed, speed up the calculation by disabling packing with pack=False.

  • 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.

  • output_prefix (str) – Prefix for output files.

results
run(start: int | None = None, stop: int | None = None, step: int | None = None, frames: int | None = None, verbose: bool | None = None, progressbar_kwargs: dict | None = None) Self

Iterate over the trajectory.

Parameters:
  • start (int) – start frame of analysis

  • stop (int) – stop frame of analysis

  • step (int) – number of frames to skip between each analysed frame

  • frames (array_like) – array of integers or booleans to slice trajectory; frames can only be used instead of start, stop, and step. Setting both frames and at least one of start, stop, step to a non-default value will raise a ValueError.

  • verbose (bool) – Turn on verbosity

  • progressbar_kwargs (dict) – ProgressBar keywords with custom parameters regarding progress bar position, etc; see MDAnalysis.lib.log.ProgressBar for full list.

Returns:

self – analysis object

Return type:

object

save() None[source]

Save results of analysis to files specified by output_prefix.