Sphere classes#
Base class for spherical analysis.
- class maicos.core.sphere.ProfileSphereBase(weighting_function: Callable, normalization: str, atomgroups: AtomGroup | List[AtomGroup], grouping: str, bin_method: str, output: str, f_kwargs: Dict | None = None, **kwargs)[source]#
Bases:
SphereBase
,ProfileBase
Base class for computing radial profiles in a spherical geometry.
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
maicos.core.base.AnalysisBase
or the General design section.- Parameters:
weighting_function (callable) – The function calculating the array weights for the histogram analysis. It must take an Atomgroup as first argument and a grouping (‘atoms’, ‘residues’, ‘segments’, ‘molecules’, ‘fragments’) as second. Additional parameters can be given as f_kwargs. The function must return a numpy.ndarray with the same length as the number of group members.
normalization ({
"none"
,"number"
,"volume"
}) – The normalization of the profile performed in every frame. IfNone
, no normalization is performed. If number, the histogram is divided by the number of occurences in each bin. If volume, the profile is divided by the volume of each bin.f_kwargs (dict) – Additional parameters for function
atomgroups (MDAnalysis.core.groups.AtomGroup or list[MDAnalysis.core.groups.AtomGroup]) – A
AtomGroup
or list thereof 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.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 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
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 Å).
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"
).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.
- results.bin_pos#
Bin positions (in Å) ranging from
rmin
tormax
.- Type:
- results.profile#
Calculated profile.
- Type:
- results.dprofile#
Estimated profile’s uncertainity.
- Type:
- class maicos.core.sphere.SphereBase(atomgroups: AtomGroup | List[AtomGroup], rmin: float, rmax: None | float, bin_width: float, **kwargs)[source]#
Bases:
AnalysisBase
Analysis class providing options and attributes for a spherical system.
Provide the results attribute r.
- Parameters:
atomgroups (MDAnalysis.core.groups.AtomGroup or list[MDAnalysis.core.groups.AtomGroup]) – A
AtomGroup
or list thereof 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.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 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
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 Å).
kwargs (dict) – Parameters parsed to AnalysisBase.
- 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)