DielectricPlanar

class maicos.DielectricPlanar(atomgroup: AtomGroup, temperature: float = 300, vcutwidth: float = 0.1, is_3d: bool = False, dim: int = 2, zmin: float | None = None, zmax: float | None = None, bin_width: float = 0.5, sym: bool = False, refgroup: AtomGroup | None = None, unwrap: bool = True, pack: bool = True, jitter: float = 0.0, concfreq: int = 0, output_prefix: str = 'eps')[source]

Bases: PlanarBase

Planar dielectric profiles.

Computes the parallel \(\varepsilon_\parallel(z)\) and inverse perpendicular (\(\varepsilon_\perp^{-1}(r)\)) components of the planar dielectric tensor \(\varepsilon\). The components are binned along the cartesian \(z\) direction yielding the component normal to the surface and defined by the dim parameter.

For usage please refer to How-to: Dielectric constant and for details on the theory see Dielectric constant measurement.

For correlation analysis, the norm of the parallel total dipole moment is used. For further information on the correlation analysis please refer to AnalysisBase or the General design section.

Also, please read and cite Schlaich et al.[1] and Refs. [2], [3].

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

  • temperature (float) – Reference temperature (K)

  • vcutwidth (float) – Spacing of virtual cuts (bins) along the parallel directions.

  • is_3d (bool) – Use 3d-periodic boundary conditions, i.e., include the dipole correction for the interaction between periodic images [4].

  • dim ({0, 1, 2}) – Dimension for binning (x=0, y=1, z=1).

  • zmin (float) – Minimal coordinate for evaluation (in Å) with respect to the center of mass of the refgroup. If zmin=None, all coordinates down to the lower cell boundary are taken into account.

  • zmax (float) – Maximal coordinate for evaluation (in Å) with respect to the center of mass of the refgroup. If zmax = None, all coordinates up to the upper cell boundary are taken into account.

  • bin_width (float) – Width of the bins (in Å).

  • sym (bool) – Symmetrize the profile. Only works in combination with refgroup.

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

Bin positions (in Å) ranging from zmin to zmax.

Type:

numpy.ndarray

results.eps_par

Reduced parallel dielectric profile \((\varepsilon_\parallel(z) - 1)\) of the selected AtomGroup

Type:

numpy.ndarray

results.deps_par

Uncertainty of parallel dielectric profile

Type:

numpy.ndarray

results.eps_par_self

Reduced self contribution of parallel dielectric profile \((\varepsilon_{\parallel,\mathrm{self}}(z) - 1)\)

Type:

numpy.ndarray

results.eps_par_coll

Reduced collective contribution of parallel dielectric profile \((\varepsilon_{\parallel,\mathrm{coll}}(z) - 1)\)

Type:

numpy.ndarray

results.eps_perp

Reduced inverse perpendicular dielectric profile \((\varepsilon^{-1}_\perp(z) - 1)\)

Type:

numpy.ndarray

results.deps_perp

Uncertainty of inverse perpendicular dielectric profile

Type:

numpy.ndarray

results.eps_perp_self

Reduced self contribution of the inverse perpendicular dielectric profile \((\varepsilon^{-1}_{\perp,\mathrm{self}}(z) - 1)\)

Type:

numpy.ndarray

results.eps_perp_coll

Reduced collective contribution of the inverse perpendicular dielectric profile \((\varepsilon^{-1}_{\perp,\mathrm{coll}}(z) - 1)\)

Type:

numpy.ndarray

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.