DensityPlanar¶
- class maicos.DensityPlanar(atomgroup: AtomGroup, dens: str = 'mass', dim: int = 2, zmin: float | None = None, zmax: float | None = None, bin_width: float = 1, bin_method: str = 'com', grouping: str = 'atoms', sym: bool = False, refgroup: AtomGroup | None = None, unwrap: bool = True, pack: bool = True, jitter: float = 0.0, concfreq: int = 0, output: str = 'density.dat')[source]¶
Bases:
ProfilePlanarBase
Cartesian partial density profiles.
Calculations are carried out for
mass
\((\rm u \cdot Å^{-3})\),number
\((\rm Å^{-3})\), partialcharge
\((\rm e \cdot Å^{-3})\) or electron \((\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 withunwrap=False
to gain a speedup. For grouping with respect to atoms, theunwrap
option is always ignored.For the correlation analysis the central bin (\(N / 2\)) 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.dens ({
"mass"
,"number"
,"charge"
,"electron"
}) – density type to be calculated.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 Å).
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 wheregrouping="residues"
,"segments"
,"molecules"
or"fragments"
).sym (bool) – Symmetrize the profile. Only works in combination with
refgroup
.refgroup (MDAnalysis.core.groups.AtomGroup) – Reference
AtomGroup
used for the calculation. Ifrefgroup
is provided, the calculation is performed relative to the center of mass of the AtomGroup. Ifrefgroup
isNone
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 useunwrap=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 everyconcfreq
frames.output (str) – Output filename.
- results.bin_pos¶
Bin positions (in Å) ranging from
zmin
tozmax
.- Type:
- results.profile¶
Calculated profile.
- Type:
- results.dprofile¶
Estimated profile’s uncertainity.
- Type:
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.
- 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 ofstart
,stop
, andstep
. Setting bothframes
and at least one ofstart
,stop
,step
to a non-default value will raise aValueError
.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: