DiporderSphere

class maicos.DiporderSphere(atomgroup: AtomGroup, order_parameter: str = 'P0', rmin: float = 0, rmax: float | None = None, bin_width: float = 1, bin_method: str = 'com', grouping: str = 'residues', refgroup: AtomGroup | None = None, unwrap: bool = True, pack: bool = True, jitter: float = 0.0, concfreq: int = 0, output: str = 'diporder_sphere.dat')[source]

Bases: ProfileSphereBase

Spherical dipolar order parameters.

Calculations include the projected dipole density \(P_0⋅ρ(z)⋅\cos(θ[z])\), the dipole orientation \(\cos(θ[z])\), the squared dipole orientation \(\cos²(Θ[z])\) and the number density \(ρ(z)\).

For the correlation analysis the 0th bin of the 0th’s group profile is used. For further information on the correlation analysis please refer to AnalysisBase or the General design section.

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

  • order_parameter ({"P0", "cos_theta", "cos_2_theta"}) –

    Order parameter to be calculated:
    • "P0": total dipole moment projected on an axis

    • "cos_theta": cosine of the dipole moment with an axis

    • "cos_2_theta": squred cosine with an axis.

  • rmin (float) – Minimal radial coordinate relative to the center of mass of the refgroup for evaluation (in Å).

  • rmax (float) –

    Maximal radial coordinate relative to the center of mass of the refgroup for evaluation (in Å).

    If rmax=None, the box extension is taken.

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

  • bin_method ({"com", "cog", "coc"}) –

    Method for the position binning.

    The possible options are center of mass ("com"), center of geometry ("cog"), and center of charge ("coc").

  • grouping ({"atoms", "residues", "segments", "molecules", "fragments"}) –

    Atom grouping for the calculations.

    The possible grouping options are the atom positions (in the case where grouping="atoms") or the center of mass of the specified grouping unit (in the case where grouping="residues", "segments", "molecules" or "fragments").

  • 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 (str) – Output filename.

results.bin_pos

Bin positions (in Å) ranging from rmin to rmax.

Type:

numpy.ndarray

results.profile

Calculated profile.

Type:

numpy.ndarray

results.dprofile

Estimated profile’s uncertainity.

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

Save results of analysis to file specified by output.