AnalysisBase¶
- class maicos.core.AnalysisBase(atomgroup: AtomGroup, unwrap: bool, pack: bool, refgroup: None | AtomGroup, jitter: float, concfreq: int, wrap_compound: str)[source]¶
Bases:
_Runner
,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.During the analysis, the correlation time of an observable can be estimated to ensure that calculated errors are reasonable. For this, the
_single_frame()
method has to return a singlefloat
. For details on the computation of the correlation and its further analysis refer tomaicos.lib.util.correlation_analysis()
.- 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.wrap_compound (str) – The group which will be kept together through the wrap processes. Allowed values are:
"atoms"
,"group"
,"residues"
,"segments"
"molecules"
"fragments". (or)
- _universe¶
The Universe the AtomGroup belong to
- _trajectory¶
The trajectory the AtomGroup 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()
- Type:
MDAnalysis.analysis.base.Results
- _obs¶
Observables of the current frame
- Type:
MDAnalysis.analysis.base.Results
- _obs.box_center¶
Center of the simulation cell of the current frame
- Type:
- sums¶
Sum of the observables across frames. Keys are the same as
_obs
.- Type:
MDAnalysis.analysis.base.Results
- sems¶
Standard errors of the mean of the observables. Keys are the same as
_obs
- Type:
MDAnalysis.analysis.base.Results
- corrtime¶
The correlation time of the analysed data. For details on how this is calculated see
maicos.lib.util.correlation_analysis()
.- Type:
- Raises:
ValueError – If any of the provided AtomGroups (atomgroup or refgroup) does not contain any atoms.
Example
To write your own analysis module you can use the example given below. As with all MAICoS modules, this inherits from the
AnalysisBase
class.The example will calculate the average box volume and stores the result within the
result
object of the class.>>> import logging >>> from typing import Optional
>>> import MDAnalysis as mda >>> import numpy as np
>>> from maicos.core import AnalysisBase >>> from maicos.lib.util import render_docs
Creating a logger makes debugging easier.
>>> logger = logging.getLogger(__name__)
In the following the analysis module itself. Due to the similar structure of all MAICoS modules you can render the parameters using the
maicos.lib.util.render_docs()
decorator. The decorator will replace special keywords with a leading$
with the actual docstring as defined inmaicos.lib.util.DOC_DICT
.>>> @render_docs ... class NewAnalysis(AnalysisBase): ... '''Analysis class calcuting the average box volume.''' ... ... def __init__( ... self, ... atomgroup: mda.AtomGroup, ... concfreq: int = 0, ... temperature: float = 300, ... output: str = "outfile.dat", ... ): ... super().__init__( ... atomgroup=atomgroup, ... refgroup=None, ... unwrap=False, ... pack=True, ... jitter=0.0, ... wrap_compound="atoms", ... concfreq=concfreq, ... ) ... ... self.temperature = temperature ... self.output = output ... ... def _prepare(self): ... '''Set things up before the analysis loop begins.''' ... # self.atomgroup refers to the provided `atomgroup` ... # self._universe refers to full universe of given `atomgroup` ... self.volume = 0 ... ... def _single_frame(self): ... '''Calculate data from a single frame of trajectory. ... ... Don't worry about normalising, just deal with a single frame. ... ''' ... # Current frame index: self._frame_index ... # Current timestep object: self._ts ... ... volume = self._ts.volume ... self.volume += volume ... ... # Eeach module should return a characteristic scalar which is used ... # by MAICoS to estimate correlations of an Analysis. ... return volume ... ... def _conclude(self): ... '''Finalise the results you've gathered. ... ... Called at the end of the run() method to finish everything up. ... ''' ... self.results.volume = self.volume / self.n_frames ... logger.info( ... "Average volume of the simulation box " ... f"{self.results.volume:.2f} ų" ... ) ... ... def save(self) -> None: ... '''Save results of analysis to file specified by ``output``. ... ... Called at the end of the run() method after _conclude. ... ''' ... self.savetxt( ... self.output, np.array([self.results.volume]), columns="volume / ų" ... ) ...
Afterwards the new analysis can be run like this
>>> import MDAnalysis as mda >>> from MDAnalysisTests.datafiles import TPR, XTC
>>> u = mda.Universe(TPR, XTC)
>>> na = NewAnalysis(u.atoms) >>> _ = na.run(start=0, stop=10) >>> print(round(na.results.volume, 2)) 362631.65
Results can also be accessed by key
>>> print(round(na.results["volume"], 2)) 362631.65
- 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 [source]¶
Iterate over the trajectory.
- savetxt(fname: str, X: ndarray, columns: List[str] | None = None) 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
atomgroup that was analyzed
output messages from modules and base classes (if they exist)