Gromacs  2020.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members | Public Member Functions | Protected Member Functions
gmx::TrajectoryAnalysisModule Class Referenceabstract

#include <gromacs/trajectoryanalysis/analysismodule.h>

Inherited by AnalysisTemplate, and gmx::SelectionTester.

Description

Base class for trajectory analysis modules.

Trajectory analysis methods should derive from this class and override the necessary virtual methods to implement initialization (initOptions(), optionsFinished(), initAnalysis(), initAfterFirstFrame()), per-frame analysis (analyzeFrame()), and final processing (finishFrames(), finishAnalysis(), writeOutput()).

For parallel analysis using threads, only a single object is constructed, but the methods startFrames(), analyzeFrame() and finishFrames() are called in each thread. Frame-local data should be initialized in startFrames() and stored in a class derived from TrajectoryAnalysisModuleData that is passed to the other methods. The default implementation of startFrames() can be used if only data handles and selections need to be thread-local.

To get the full benefit from this class, analysis data objects and selections should be used in the implementation. See the corresponding modules' documentation for details of how they work.

Typical way of using AnalysisData in derived classes is to have the AnalysisData object as a member variable and register it using registerAnalysisDataset(). Analysis modules are initialized in initAnalysis() and the processing chain is initialized. If any of the modules is required, e.g., for post-processing in finishAnalysis(), it can be stored in a member variable. To add data to the data object in analyzeFrame(), a data handle is obtained using TrajectoryAnalysisModuleData::dataHandle().

Typical way of using selections in derived classes is to have the required Selection objects (or SelectionList objects) as member variables, and add the required selection options in initOptions(). These member variables can be accessed in initAnalysis() to get general information about the selections. In analyzeFrame(), these selection objects should not be used directly, but instead TrajectoryAnalysisModuleData::parallelSelection() should be used to obtain a selection object that works correctly also for parallel analysis.

Derived classes should use exceptions to indicate errors in the virtual methods.

Examples:
template.cpp.

Public Member Functions

virtual void initOptions (IOptionsContainer *options, TrajectoryAnalysisSettings *settings)=0
 Initializes options understood by the module. More...
 
virtual void optionsFinished (TrajectoryAnalysisSettings *settings)
 Called after all option values have been set. More...
 
virtual void initAnalysis (const TrajectoryAnalysisSettings &settings, const TopologyInformation &top)=0
 Initializes the analysis. More...
 
virtual void initAfterFirstFrame (const TrajectoryAnalysisSettings &settings, const t_trxframe &fr)
 Performs additional initialization after reading the first frame. More...
 
virtual
TrajectoryAnalysisModuleDataPointer 
startFrames (const AnalysisDataParallelOptions &opt, const SelectionCollection &selections)
 Starts the analysis of frames. More...
 
virtual void analyzeFrame (int frnr, const t_trxframe &fr, t_pbc *pbc, TrajectoryAnalysisModuleData *pdata)=0
 Analyzes a single frame. More...
 
virtual void finishFrames (TrajectoryAnalysisModuleData *pdata)
 Finishes the analysis of frames. More...
 
virtual void finishAnalysis (int nframes)=0
 Postprocesses data after frames have been read. More...
 
virtual void writeOutput ()=0
 Writes output into files and/or standard output/error. More...
 
int datasetCount () const
 Returns the number of datasets provided by the module. More...
 
const std::vector< std::string > & datasetNames () const
 Returns a vector with the names of datasets provided by the module. More...
 
AbstractAnalysisDatadatasetFromIndex (int index) const
 Returns a pointer to the data set index. More...
 
AbstractAnalysisDatadatasetFromName (const char *name) const
 Returns a pointer to the data set with name name. More...
 
void finishFrameSerial (int frameIndex)
 Processes data in AnalysisData objects in serial for each frame. More...
 

Protected Member Functions

 TrajectoryAnalysisModule ()
 Initializes the dataset registration mechanism. More...
 
void registerBasicDataset (AbstractAnalysisData *data, const char *name)
 Registers a dataset that exports data. More...
 
void registerAnalysisDataset (AnalysisData *data, const char *name)
 Registers a parallelized dataset that exports data. More...
 

Constructor & Destructor Documentation

gmx::TrajectoryAnalysisModule::TrajectoryAnalysisModule ( )
protected

Initializes the dataset registration mechanism.

Exceptions
std::bad_allocif out of memory.

Member Function Documentation

virtual void gmx::TrajectoryAnalysisModule::analyzeFrame ( int  frnr,
const t_trxframe &  fr,
t_pbc pbc,
TrajectoryAnalysisModuleData pdata 
)
pure virtual

Analyzes a single frame.

Parameters
[in]frnrFrame number, a zero-based index that uniquely identifies the frame.
[in]frCurrent frame.
[in]pbcPeriodic boundary conditions for fr.
[in,out]pdataData structure for frame-local data.

This method is called once for each frame to be analyzed, and should analyze the positions provided in the selections. Data handles and selections should be obtained from the pdata structure.

For threaded analysis, this method is called asynchronously in different threads to analyze different frames. The pdata structure is one of the structures created with startFrames(), but no assumptions should be made about which of these data structures is used. It is guaranteed that two instances of analyzeFrame() are not running concurrently with the same pdata data structure. Any access to data structures not stored in pdata should be designed to be thread-safe.

int gmx::TrajectoryAnalysisModule::datasetCount ( ) const

Returns the number of datasets provided by the module.

Does not throw.

AbstractAnalysisData & gmx::TrajectoryAnalysisModule::datasetFromIndex ( int  index) const

Returns a pointer to the data set index.

Parameters
[in]indexData set to query for.
Returns
Reference to the requested data set.
Exceptions
APIErrorif index is not valid.

index should be >= 0 and < datasetCount().

The return value is not const to allow callers to add modules to the data sets. However, the AbstractAnalysisData interface does not provide any means to alter the data, so the module does not need to care about external modifications.

AbstractAnalysisData & gmx::TrajectoryAnalysisModule::datasetFromName ( const char *  name) const

Returns a pointer to the data set with name name.

Parameters
[in]nameData set to query for.
Returns
Reference to the requested data set.
Exceptions
APIErrorif name is not valid.

name should be one of the names returned by datasetNames().

The return value is not const to allow callers to add modules to the data sets. However, the AbstractAnalysisData interface does not provide any means to alter the data, so the module does not need to care about external modifications.

const std::vector< std::string > & gmx::TrajectoryAnalysisModule::datasetNames ( ) const

Returns a vector with the names of datasets provided by the module.

Does not throw.

virtual void gmx::TrajectoryAnalysisModule::finishAnalysis ( int  nframes)
pure virtual

Postprocesses data after frames have been read.

Parameters
[in]nframesTotal number of frames processed.

This function is called after all finishFrames() calls have been called. nframes will equal the number of calls to analyzeFrame() that have occurred.

void gmx::TrajectoryAnalysisModule::finishFrames ( TrajectoryAnalysisModuleData pdata)
virtual

Finishes the analysis of frames.

Parameters
[in]pdataData structure for thread-local data.

This method is called once for each call of startFrames(), with the data structure returned by the corresponding startFrames(). The pdata object should be destroyed by the caller after this function has been called.

You only need to override this method if you need custom operations to combine data from the frame-local data structures to get the final result. In such cases, the data should be aggregated in this function and stored in a member attribute.

The default implementation does nothing.

See Also
startFrames()
void gmx::TrajectoryAnalysisModule::finishFrameSerial ( int  frameIndex)

Processes data in AnalysisData objects in serial for each frame.

Parameters
[in]frameIndexIndex of the frame that has been finished.

This method is called by the framework in order for each frame, after the analysis for that frame has been finished. These calls always execute in serial and in sequential frame order, even during parallel analysis where multiple analyzeFrame() calls may be executing concurrently.

See Also
AnalysisData::finishFrameSerial()
void gmx::TrajectoryAnalysisModule::initAfterFirstFrame ( const TrajectoryAnalysisSettings settings,
const t_trxframe &  fr 
)
virtual

Performs additional initialization after reading the first frame.

When this function is called, selections are the same as in initAnalysis(), i.e., they have not been evaluated for the first frame.

It is necessary to override this method only if the module needs to do initialization for which it requires data from the first frame.

The default implementation does nothing.

virtual void gmx::TrajectoryAnalysisModule::initAnalysis ( const TrajectoryAnalysisSettings settings,
const TopologyInformation &  top 
)
pure virtual

Initializes the analysis.

Parameters
[in]settingsSettings to pass to and from the module.
[in]topTopology information.

When this function is called, selections have been initialized based on user input, and a topology has been loaded if provided by the user. For dynamic selections, the selections have been evaluated to the largest possible selection, i.e., the selections passed to analyzeFrame() are always a subset of the selections provided here.

virtual void gmx::TrajectoryAnalysisModule::initOptions ( IOptionsContainer options,
TrajectoryAnalysisSettings settings 
)
pure virtual

Initializes options understood by the module.

Parameters
[in,out]optionsOptions object to add the options to.
[in,out]settingsSettings to pass to and from the module.

This method is called first after the constructor, and it should add options understood by the module to options. Output values from options (including selections) should be stored in member variables.

In addition to initializing the options, this method can also provide information about the module's requirements using the settings object; see TrajectoryAnalysisSettings for more details.

If settings depend on the option values provided by the user, see optionsFinished().

void gmx::TrajectoryAnalysisModule::optionsFinished ( TrajectoryAnalysisSettings settings)
virtual

Called after all option values have been set.

Parameters
[in,out]settingsSettings to pass to and from the module.

This method is called after option values have been assigned (but interactive selection input has not yet been performed).

If the module needs to change settings that affect topology loading (can be done using the settings object) or selection initialization (can be done using SelectionOptionInfo) based on option values, this method has to be overridden.

The default implementation does nothing.

void gmx::TrajectoryAnalysisModule::registerAnalysisDataset ( AnalysisData data,
const char *  name 
)
protected

Registers a parallelized dataset that exports data.

Parameters
dataAnalysisData object to register.
[in]nameName to register the dataset with.
Exceptions
std::bad_allocif out of memory.

This method works as registerBasicDataset(), but additionally allows data handles for data to be accessed using TrajectoryAnalysisData.

See Also
registerBasicDataset()
void gmx::TrajectoryAnalysisModule::registerBasicDataset ( AbstractAnalysisData data,
const char *  name 
)
protected

Registers a dataset that exports data.

Parameters
dataData object to register.
[in]nameName to register the dataset with.
Exceptions
std::bad_allocif out of memory.

Registers data as a dataset that provides output from the analysis module. Callers for the module can access the dataset with datasetFromName() using name as an AbstractAnalysisData object. This allows them to add their own data modules to do extra processing.

name must be unique across all calls within the same TrajectoryAnalysisModule instance.

TrajectoryAnalysisModuleDataPointer gmx::TrajectoryAnalysisModule::startFrames ( const AnalysisDataParallelOptions &  opt,
const SelectionCollection selections 
)
virtual

Starts the analysis of frames.

Parameters
[in]optParallel options
[in]selectionsFrame-local selection collection object.
Returns
Data structure for thread-local data.

This function is necessary only for threaded parallelization. It is called once for each thread and should initialize a class that contains any required frame-local data in the returned value. The default implementation creates a basic data structure that holds thread-local data handles for all data objects registered with registerAnalysisDataset(), as well as the thread-local selection collection. These can be accessed in analyzeFrame() using the methods in TrajectoryAnalysisModuleData. If other thread-local data is needed, this function should be overridden and it should create an instance of a class derived from TrajectoryAnalysisModuleData.

See Also
TrajectoryAnalysisModuleData
virtual void gmx::TrajectoryAnalysisModule::writeOutput ( )
pure virtual

Writes output into files and/or standard output/error.

All output from the module, excluding data written out for each frame during analyzeFrame(), should be confined into this function. This function is guaranteed to be called only after finishAnalysis().


The documentation for this class was generated from the following files: