Merge "Improve selection interface (remove some pointers)."
[alexxy/gromacs.git] / src / gromacs / trajectoryanalysis / analysismodule.h
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11  * Copyright (c) 2001-2009, The GROMACS development team,
12  * check out http://www.gromacs.org for more information.
13
14  * This program is free software; you can redistribute it and/or
15  * modify it under the terms of the GNU General Public License
16  * as published by the Free Software Foundation; either version 2
17  * of the License, or (at your option) any later version.
18  *
19  * If you want to redistribute modifications, please consider that
20  * scientific software is very special. Version control is crucial -
21  * bugs must be traceable. We will be happy to consider code for
22  * inclusion in the official distribution, but derived work must not
23  * be called official GROMACS. Details are found in the README & COPYING
24  * files - if they are missing, get the official version at www.gromacs.org.
25  *
26  * To help us fund GROMACS development, we humbly ask that you cite
27  * the papers on the package - you can find them in the top README file.
28  *
29  * For more info, check our website at http://www.gromacs.org
30  */
31 /*! \file
32  * \brief
33  * Declares gmx::TrajectoryAnalysisModule and
34  * gmx::TrajectoryAnalysisModuleData.
35  *
36  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
37  * \inpublicapi
38  * \ingroup module_trajectoryanalysis
39  */
40 #ifndef GMX_TRAJECTORYANALYSIS_ANALYSISMODULE_H
41 #define GMX_TRAJECTORYANALYSIS_ANALYSISMODULE_H
42
43 #include <string>
44 #include <vector>
45
46 #include "../legacyheaders/typedefs.h"
47
48 #include "../selection/selection.h" // For gmx::SelectionList
49 #include "../utility/common.h"
50 #include "../utility/uniqueptr.h"
51
52 namespace gmx
53 {
54
55 class AbstractAnalysisData;
56 class AnalysisData;
57 class AnalysisDataHandle;
58 class AnalysisDataParallelOptions;
59 class Options;
60 class SelectionCollection;
61 class TopologyInformation;
62 class TrajectoryAnalysisModule;
63 class TrajectoryAnalysisSettings;
64
65 /*! \brief
66  * Base class for thread-local data storage during trajectory analysis.
67  *
68  * Thread-local storage of data handles and selections is implemented in this
69  * class; TrajectoryAnalysisModule instances can access the thread-local values
70  * using dataHandle() and parallelSelection().
71  *
72  * \see TrajectoryAnalysisModule::startFrames()
73  * \see TrajectoryAnalysisModule::analyzeFrame()
74  * \see TrajectoryAnalysisModule::finishFrames()
75  *
76  * \inpublicapi
77  * \ingroup module_trajectoryanalysis
78  */
79 class TrajectoryAnalysisModuleData
80 {
81     public:
82         virtual ~TrajectoryAnalysisModuleData();
83
84         /*! \brief
85          * Performs any finishing actions after all frames have been processed.
86          *
87          * This function is called immediately before the destructor.
88          * All implementations should call finishDataHandles().
89          */
90         virtual void finish() = 0;
91
92         /*! \brief
93          * Returns a data handle for a given dataset.
94          *
95          * Allowed data sets are those that have been registered with
96          * TrajectoryAnalysisModule::registerAnalysisDataset().
97          */
98         AnalysisDataHandle dataHandle(const AnalysisData &data);
99         /*! \brief
100          * Returns a selection that corresponds to the given selection.
101          */
102         Selection parallelSelection(const Selection &selection);
103         /*! \brief
104          * Returns a set of selection that corresponds to the given selections.
105          */
106         SelectionList parallelSelections(const SelectionList &selections);
107
108     protected:
109         /*! \brief
110          * Initializes thread-local storage for data handles and selections.
111          *
112          * \param[in] module     Analysis module to use for data objects.
113          * \param[in] opt        Data parallelization options.
114          * \param[in] selections Thread-local selection collection.
115          *
116          * Calls AnalysisData::startData() on all data objects registered with
117          * TrajectoryAnalysisModule::registerAnalysisDataset() in \p module.
118          * The handles are accessible through dataHandle().
119          */
120         TrajectoryAnalysisModuleData(TrajectoryAnalysisModule *module,
121                                      const AnalysisDataParallelOptions &opt,
122                                      const SelectionCollection &selections);
123
124         /*! \brief
125          * Calls finishData() on all data handles.
126          *
127          * This function should be called from the implementation of finish()
128          * in all subclasses.
129          */
130         void finishDataHandles();
131
132     private:
133         class Impl;
134
135         PrivateImplPointer<Impl> _impl;
136 };
137
138 //! Smart pointer to manage a TrajectoryAnalysisModuleData object.
139 typedef gmx_unique_ptr<TrajectoryAnalysisModuleData>::type
140         TrajectoryAnalysisModuleDataPointer;
141
142 /*! \brief
143  * Base class for trajectory analysis methods.
144  *
145  * Trajectory analysis methods should derive from this class and override the
146  * necessary virtual functions to implement initialization (initOptions(),
147  * initOptionsDone(), initAnalysis(), initAfterFirstFrame()), per-frame analysis
148  * (analyzeFrame()), and final processing (finishFrames(), finishAnalysis(),
149  * writeOutput()).
150  *
151  * For parallel analysis using threads, only a single object is constructed,
152  * but the methods startFrames(), analyzeFrame() and finishFrames() are called
153  * in each thread.  Frame-local data should be initialized in startFrames() and
154  * stored in a class derived from TrajectoryAnalysisModuleData that is passed
155  * to the other methods.  The default implementation of startFrames() can be
156  * used if only data handles and selections need to be thread-local.
157  *
158  * \inpublicapi
159  * \ingroup module_trajectoryanalysis
160  */
161 class TrajectoryAnalysisModule
162 {
163     public:
164         virtual ~TrajectoryAnalysisModule();
165
166         /*! \brief
167          * Initializes options understood by the module.
168          *
169          * In addition to initializing the options, this function can also
170          * provide information about its requirements using the \p settings
171          * object; see TrajectoryAnalysisSettings for more details.
172          *
173          * If settings depend on the option values provided by the user, see
174          * initOptionsDone().
175          */
176         virtual Options &initOptions(TrajectoryAnalysisSettings *settings) = 0;
177         /*! \brief
178          * Called after all option values have been set.
179          *
180          * If the module needs to change settings that affect topology loading
181          * or selection initialization based on option values, this function
182          * has to be overridden.
183          *
184          * The default implementation does nothing.
185          */
186         virtual void initOptionsDone(TrajectoryAnalysisSettings *settings);
187         /*! \brief
188          * Initializes the analysis.
189          *
190          * When this function is called, selections have been initialized based
191          * on user input, and a topology has been loaded if provided by the
192          * user.  For dynamic selections, the selections have been evaluated to
193          * the largest possible selection, i.e., the selections passed to
194          * analyzeFrame() are always a subset of the selections provided here.
195          */
196         virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
197                                   const TopologyInformation &top) = 0;
198         /*! \brief
199          * Performs additional initialization after reading the first frame.
200          *
201          * When this function is called, selections are the same as in
202          * initAnalysis(), i.e., they have not been evaluated for the first
203          * frame.
204          *
205          * It is necessary to override this method only if the module needs to
206          * do initialization for which it requires data from the first frame.
207          *
208          * The default implementation does nothing.
209          */
210         virtual void initAfterFirstFrame(const t_trxframe &fr);
211
212         /*! \brief
213          * Starts the analysis of frames.
214          *
215          * \param[in]  opt
216          * \param[in]  selections  Frame-local selection collection object.
217          * \returns    Data structure for thread-local data.
218          *
219          * This function is necessary only for threaded parallelization.
220          * It is called once for each thread and should initialize a class that
221          * contains any required frame-local data in the returned value.
222          * The default implementation creates a basic data structure that holds
223          * thread-local data handles for all data objects registered with
224          * registerAnalysisDataset(), as well as the thread-local selection
225          * collection.  These can be accessed in analyzeFrame() using the
226          * methods in TrajectoryAnalysisModuleData.
227          * If other thread-local data is needed, this function should be
228          * overridden and it should create an instance of a class derived from
229          * TrajectoryAnalysisModuleData.
230          *
231          * \see TrajectoryAnalysisModuleData
232          */
233         virtual TrajectoryAnalysisModuleDataPointer startFrames(
234                 const AnalysisDataParallelOptions &opt,
235                 const SelectionCollection &selections);
236         /*! \brief
237          * Analyzes a single frame.
238          *
239          * \param[in]     frnr   Frame number, a zero-based index that
240          *      uniquely identifies the frame.
241          * \param[in]     fr     Current frame.
242          * \param[in]     pbc    Periodic boundary conditions for \p fr.
243          * \param[in,out] pdata  Data structure for frame-local data.
244          *
245          * This function is called once for each frame to be analyzed,
246          * and should analyze the positions provided in \p sel.
247          *
248          * For threaded analysis, this function is called asynchronously in
249          * different threads to analyze different frames. The \p pdata
250          * structure is one of the structures created with startFrames(),
251          * but no assumptions should be made about which of these data
252          * structures is used. It is guaranteed that two instances of
253          * analyzeFrame() are not running concurrently with the same \p pdata
254          * data structure.
255          * Any access to data structures not stored in \p pdata should be
256          * designed to be thread-safe.
257          */
258         virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
259                                   TrajectoryAnalysisModuleData *pdata) = 0;
260         /*! \brief
261          * Finishes the analysis of frames.
262          *
263          * \param[in]  pdata    Data structure for thread-local data.
264          *
265          * This function is called once for each call of startFrames(),
266          * with the data structure returned by the corresponding startFrames().
267          * The \p pdata object should be destroyed by the caller after this
268          * function has been called.
269          *
270          * You only need to override this method if you need custom
271          * operations to combine data from the frame-local data structures
272          * to get the final result. In such cases, the data should be
273          * aggregated in this function and stored in a member attribute.
274          *
275          * The default implementation does nothing.
276          *
277          * \see startFrames()
278          */
279         virtual void finishFrames(TrajectoryAnalysisModuleData *pdata);
280
281         /*! \brief
282          * Postprocesses data after frames have been read.
283          *
284          * This function is called after all finishFrames() calls have been
285          * called.
286          */
287         virtual void finishAnalysis(int nframes) = 0;
288         /*! \brief
289          * Writes output into files and/or standard output/error.
290          *
291          * All output from the module, excluding data written out for each
292          * frame during analyzeFrame(), should be confined into this function.
293          * This function is guaranteed to be called only after
294          * finishAnalysis().
295          */
296         virtual void writeOutput() = 0;
297
298         /*! \brief
299          * Returns the number of datasets provided by the module.
300          */
301         int datasetCount() const;
302         /*! \brief
303          * Returns a vector with the names of the datasets.
304          */
305         const std::vector<std::string> &datasetNames() const;
306         /*! \brief
307          * Returns a pointer to the data set \p index.
308          *
309          * \param[in] index  Data set to query for.
310          * \returns   A pointer to the data set, or NULL if \p index is not
311          *      valid.
312          *
313          * The return value is not const to allow callers to add modules to the
314          * data sets. However, the AbstractAnalysisData interface does not
315          * provide any means to alter the data, so the module does not need to
316          * care about external modifications.
317          */
318         AbstractAnalysisData &datasetFromIndex(int index) const;
319         /*! \brief
320          * Returns a pointer to the data set with name \p name
321          *
322          * \param[in] name  Data set to query for.
323          * \returns   A pointer to the data set, or NULL if \p name is not
324          *      recognized.
325          *
326          * The return value is not const to allow callers to add modules to the
327          * data sets. However, the AbstractAnalysisData interface does not
328          * provide any means to alter the data, so the module does not need to
329          * care about external modifications.
330          */
331         AbstractAnalysisData &datasetFromName(const char *name) const;
332
333     protected:
334         //! Initializes the dataset registration mechanism.
335         TrajectoryAnalysisModule();
336
337         /*! \brief
338          * Registers a dataset that exports data.
339          */
340         void registerBasicDataset(AbstractAnalysisData *data, const char *name);
341         /*! \brief
342          * Registers a parallelized dataset that exports data.
343          */
344         void registerAnalysisDataset(AnalysisData *data, const char *name);
345
346     private:
347         class Impl;
348
349         PrivateImplPointer<Impl> _impl;
350
351         /*! \brief
352          * Needed to access the registered analysis data sets.
353          */
354         friend class TrajectoryAnalysisModuleData;
355 };
356
357 //! Smart pointer to manage a TrajectoryAnalysisModule.
358 typedef gmx_unique_ptr<TrajectoryAnalysisModule>::type
359         TrajectoryAnalysisModulePointer;
360
361 } // namespace gmx
362
363 #endif