DielectricSpectrum Module

Module for computing dielectric spectra for bulk systems.

class spectrumkit.dielectricspectrum.Current(atomgroup: AtomGroup, unwrap: bool = False, pack: bool = False, concfreq: int = 0, jitter: float = 0.0)[source]

Bases: AnalysisBase

Current time series from md trajectory.

This module, given a molecular dynamics trajectory, produces a time series of the curren of the atomgroup and the whole system.

\(\mathbf{P} = \sum_i q_i \mathbf{r}_i\), where \(q_i\) and \(\mathbf{r}_i\) are the charge and position of atom \(i\), respectively.

Please read and cite Rinne paper…

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.

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

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

  • temperature (float) – Reference temperature (K)

  • output_prefix (str) – Prefix for output files.

results

References

save() None[source]

Save results of analysis to file specified by output.

class spectrumkit.dielectricspectrum.DielectricSpectrum(atomgroup: AtomGroup, refgroup: AtomGroup | None = None, unwrap: bool = True, pack: bool = True, concfreq: int = 0, temperature: float = 300, output_prefix: str = '', segs: int = 20, df: float | None = None, bins: int = 200, binafter: float = 20, nobin: bool = False, jitter: float = 0.0)[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.

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

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

  • temperature (float) – Reference temperature (K)

  • output_prefix (str) – Prefix for output files.

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

results

References

save() None[source]

Save results of analysis to file specified by output.

class spectrumkit.dielectricspectrum.Polarization(atomgroup: AtomGroup, unwrap: bool = True, pack: bool = False, concfreq: int = 0, jitter: float = 0.0)[source]

Bases: AnalysisBase

Dipole moment time series from md trajectory.

This module, given a molecular dynamics trajectory, produces a time series of the total dipole moment of the atomgroup and the whole system. It repairs molecules that are broken across periodic boundaries by reassembling them before calculating the dipole moment. The dipole moment is calculated as \(\mathbf{P} = \sum_i q_i \mathbf{r}_i\), where \(q_i\) and \(\mathbf{r}_i\) are the charge and position of atom \(i\), respectively.

Please read and cite [1].

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.

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

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

  • temperature (float) – Reference temperature (K)

  • output_prefix (str) – Prefix for output files.

results

References

save() None[source]

Save results of analysis to file specified by output.

class spectrumkit.dielectricspectrum.VACF(atomgroup: AtomGroup, unwrap: bool = True, pack: bool = False, concfreq: int = 0, jitter: float = 0.0)[source]

Bases: AnalysisBase

Velocity autocorrelation function from md trajectory.

This module, given a molecular dynamics trajectory, produces a time series of the velocity autocorrelation function of the atomgroup and the whole system. It repairs molecules that are broken across periodic boundaries by reassembling them before calculating the dipole moment. The dipole moment is calculated as \(\mathbf{P} = \sum_i q_i \mathbf{r}_i\), where \(q_i\) and \(\mathbf{r}_i\) are the charge and position of atom \(i\), respectively.

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.

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

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

  • temperature (float) – Reference temperature (K)

  • output_prefix (str) – Prefix for output files.

results

References

save() None[source]

Save results of analysis to file specified by output.

spectrumkit.dielectricspectrum.WienerKhinchin(time, timeseries, volume, temperature, spectrum_type='dipole')[source]

Calculate dielectric spectrum from timeseries using Wiener-Khinchin theorem.

Parameters:
  • time (np.ndarray) – Time steps of the data (n_frames) in ps.

  • timeseries (np.ndarray) – Time series data with shape (n_frames, n_coordinates).

  • volume (float) – System volume in Angstrom^3.

  • temperature (float) – System temperature in Kelvin.

  • spectrum_type (str, default="dipole") – Type of spectrum to calculate. Either “dipole” or “flux”.

Returns:

Tuple of (frequencies, complex dielectric spectrum).

Return type:

tuple[np.ndarray, np.ndarray]

Notes

Units: - dipole: e·Å - flux: e·Å / ps - dt: ps - volume: ų

References

https://dx.doi.org/10.1021/acs.jpca.0c04063

spectrumkit.dielectricspectrum.calculate_spectrum_from_current(current: ndarray, dt: float, volume: float, temperature: float, segs: int | None = None, df: float | None = None, bins: int = 200, binafter: float = 20, nobin: bool = False) dict[str, ndarray][source]

Calculate dielectric spectrum from current (dipole derivative) time series.

This function computes the complex dielectric susceptibility from a current time series J(t) = dP/dt using the Fluctuation-Dissipation theorem. It is equivalent to calculate_spectrum_from_dipole() but works directly with the time derivative of the dipole moment, avoiding the need to integrate the current back to a dipole moment.

Since FT(dP/dt) = iν·FT(P), we have \|FT(J)\|² = ν²·\|FT(P)\|², so \|FT(J)\|²/ν yields the same spectrum as \|FT(P)\|²·ν.

Parameters:
  • current (np.ndarray) – Current (dipole derivative) time series with shape (n_frames, 3) in units of e·Å/ps.

  • dt (float) – Time step between frames in picoseconds.

  • volume (float) – Average system volume in Angstrom^3.

  • temperature (float) – System temperature in Kelvin.

  • segs (int, optional) – Number of segments to break the trajectory into. If None and df is None, defaults to 20.

  • df (float, optional) – Desired frequency spacing in THz. Overrides segs if provided.

  • bins (int, default=200) – Number of bins for data averaging (logarithmic binning).

  • binafter (float, default=20) – Number of low-frequency data points left unbinned.

  • nobin (bool, default=False) – If True, prevents data binning.

Returns:

Dictionary containing: - ‘t’: time array - ‘nu’: frequency array (THz) - ‘susc’: complex susceptibility - ‘dsusc’: standard deviation of susceptibility - ‘nu_binned’: binned frequencies (if binning applied) - ‘susc_binned’: binned susceptibility (if binning applied) - ‘dsusc_binned’: binned std deviation (if binning applied)

Return type:

dict[str, np.ndarray]

spectrumkit.dielectricspectrum.calculate_spectrum_from_dipole(dipole_moment: ndarray, dt: float, volume: float, temperature: float, segs: int | None = None, df: float | None = None, bins: int = 200, binafter: float = 20, nobin: bool = False) dict[str, ndarray][source]

Calculate dielectric spectrum from dipole moment time series.

This function computes the complex dielectric susceptibility from a dipole moment time series using the Fluctuation-Dissipation theorem. It can be used independently of MDAnalysis trajectories.

Parameters:
  • dipole_moment (np.ndarray) – Dipole moment time series with shape (n_frames, 3) in units of e·Å.

  • dt (float) – Time step between frames in picoseconds.

  • volume (float) – Average system volume in Angstrom^3.

  • temperature (float) – System temperature in Kelvin.

  • segs (int, optional) – Number of segments to break the trajectory into. If None and df is None, defaults to 20.

  • df (float, optional) – Desired frequency spacing in THz. Overrides segs if provided.

  • bins (int, default=200) – Number of bins for data averaging (logarithmic binning).

  • binafter (float, default=20) – Number of low-frequency data points left unbinned.

  • nobin (bool, default=False) – If True, prevents data binning.

Returns:

Dictionary containing: - ‘t’: time array - ‘nu’: frequency array (THz) - ‘susc’: complex susceptibility - ‘dsusc’: standard deviation of susceptibility - ‘nu_binned’: binned frequencies (if binning applied) - ‘susc_binned’: binned susceptibility (if binning applied) - ‘dsusc_binned’: binned std deviation (if binning applied)

Return type:

dict[str, np.ndarray]

Notes

The algorithm is based on the Fluctuation Dissipation Relation: χ(f) = -1/(3 V k_B T ε_0) L[θ(t) ⟨P(0) dP(t)/dt⟩] where L is the Laplace transformation.