Base class#
Base class for building Analysis classes.
- class maicos.core.base.AnalysisBase(atomgroups, multi_group=False, refgroup=None, unwrap=False, jitter=None, concfreq=0)[source]#
Bases:
AnalysisBase
Base class derived from MDAnalysis for defining multi-frame analysis.
The class is designed as a template for creating multi-frame analyses. This class will automatically take care of setting up the trajectory reader for iterating, and it offers to show a progress meter. Computed results are stored inside the
results
attribute. To define a new analysis, AnalysisBase needs to be subclassed and_single_frame()
must be defined. It is also possible to define_prepare()
and_conclude()
for pre- and post-processing. All results should be stored as attributes of theMDAnalysis.analysis.base.Results
container.- Parameters:
atomgroups (list[AtomGroup]) – a list of
AtomGroup
objects for which the calculations are performed.refgroup (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.unwrap (bool) –
When
unwrap = 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.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
multi_group (bool) – Analysis is able to work with list of atomgroups
- atomgroup#
Atomgroup taken for the Analysis (available if multi_group = False)
- atomgroups#
Atomgroups taken for the Analysis (available if multi_group = True)
- _universe#
The Universe the atomgroups belong to
- _trajectory#
The trajectory the atomgroups belong to
- times#
array of Timestep times. Only exists after calling
AnalysisBase.run()
- Type:
- frames#
array of Timestep frame indices. Only exists after calling
AnalysisBase.run()
- Type:
- results#
results of calculation are stored after call to
AnalysisBase.run()
- _obs#
Observables of the current frame
- _obs.box_center#
Center of the simulation cell of the current frame
- Type:
- Raises:
ValueError – If any of the provided AtomGroups (atomgroups or refgroup) does not contain any atoms.
- property box_center#
Center of the simulation cell.
- savetxt(fname, X, columns=None)[source]#
Save to text.
An extension of the numpy savetxt function. Adds the command line input to the header and checks for a doubled defined filesuffix.
Return a header for the text file to save the data to. This method builds a generic header that can be used by any MAICoS module. It is called by the save method of each module.
- The information it collects is:
timestamp of the analysis
name of the module
version of MAICoS that was used
command line arguments that were used to run the module
module call including the default arguments
number of frames that were analyzed
atomgroups that were analyzed
output messages from modules and base classes (if they exist)
- class maicos.core.base.ProfileBase(atomgroups, weighting_function, normalization, grouping, bin_method, output, f_kwargs=None)[source]#
Bases:
object
Base class for computing profiles.
- Parameters:
atomgroups (list[AtomGroup]) – a list of
AtomGroup
objects for which the calculations are performed.grouping (str {
'atoms'
,'residues'
,'segments'
,'molecules'
,'fragments'
}) –Atom grouping for the calculations of profiles.
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 (str {
'cog'
,'com'
,'coc'
}) –Method for the position binning.
The possible options are center of geometry (
cog
), center of mass (com
), and center of charge (coc
).output (str) – Output filename.
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 (str {'None', 'number', 'volume'}) – The normalization of the profile performed in every frame. If None 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
- results.profile#
Calculated profile.
- Type:
- results.dprofile#
Estimated profile’s uncertainity.
- Type: