DensityPlanar#

class maicos.DensityPlanar(atomgroup: AtomGroup, dens: str = 'mass', dim: int = 2, zmin: float | None = None, zmax: float | None = None, bin_width: float = 1, refgroup: AtomGroup | None = None, sym: bool = False, grouping: str = 'atoms', unwrap: bool = True, bin_method: str = 'com', output: str = 'density.dat', concfreq: int = 0, jitter: float = 0.0)[source]#

Bases: ProfilePlanarBase

Cartesian partial density profiles.

Calculations are carried out for mass \((\rm u \cdot Å^{-3})\), number \((\rm Å^{-3})\) or charge \((\rm e \cdot Å^{-3})\) density profiles along certain cartesian axes [x, y, z] of the simulation cell. Cell dimensions are allowed to fluctuate in time.

For grouping with respect to molecules, residues etc., the corresponding centers (i.e., center of mass), taking into account periodic boundary conditions, are calculated. For these calculations molecules will be unwrapped/made whole. Trajectories containing already whole molecules can be run with unwrap=False to gain a speedup. For grouping with respect to atoms, the unwrap option is always ignored.

For the correlation analysis the central bin (\(N \backslash 2\)) of the 0th’s group profile is used. For further information on the correlation analysis please refer to maicos.core.base.AnalysisBase or the General design section.

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.

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

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

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

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

  • output (str) – Output filename.

  • dens ({"mass", "number", "charge"}) – density type to be calculated.

results.bin_pos#

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

Type:

numpy.ndarray

results.profile#

Calculated profile.

Type:

numpy.ndarray

results.dprofile#

Estimated profile’s uncertainity.

Type:

numpy.ndarray

Notes

Partial mass density profiles can be used to calculate the ideal component of the chemical potential. For details, take a look at the corresponding How-to guide.