SphereBase¶
- class maicos.core.SphereBase(atomgroup: AtomGroup, unwrap: bool, pack: bool, refgroup: AtomGroup | None, jitter: float, concfreq: int, rmin: float, rmax: None | float, bin_width: float, wrap_compound: str)[source]¶
Bases:
AnalysisBase
Analysis class providing options and attributes for a spherical system.
Provide the results attribute r.
- 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 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
.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.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. 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 Å).wrap_compound (str) – The group which will be kept together through the wrap processes. Allowed values are:
"atoms"
,"group"
,"residues"
,"segments"
"molecules"
"fragments". (or)
- results.bin_pos¶
Bin positions (in Å) ranging from
rmin
tormax
.- Type:
- pos_sph¶
positions in spherical coordinats (r, phi, theta)
- Type:
- _obs.bin_pos¶
Central bin position of each bin (in Å) in the current frame.
- Type:
numpy.ndarray, (n_bins)
- _obs.bin_edges¶
Edges of the bins (in Å) in the current frame.
- Type:
numpy.ndarray, (n_bins + 1)
- _obs.bin_area¶
Surface area (in Å^2) of the sphere of each bin with radius bin_pos in the current frame. Calculated via \(4 \pi r_i^2\) where \(i\) is the index of the bin.
- Type:
numpy.ndarray, (n_bins)
- results.bin_volume¶
volume of a spherical shell of each bins (in Å^3) of the current frame. Calculated via \(4\pi/3 \left(r_{i+1}^3 - r_i^3 \right)\) where i is the index of the bin.
- Type:
numpy.ndarray, (n_bins)