Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / trajectoryanalysis / analysismodule.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \file
36  * \brief
37  * Declares gmx::TrajectoryAnalysisModule and
38  * gmx::TrajectoryAnalysisModuleData.
39  *
40  * \author Teemu Murtola <teemu.murtola@gmail.com>
41  * \inpublicapi
42  * \ingroup module_trajectoryanalysis
43  */
44 #ifndef GMX_TRAJECTORYANALYSIS_ANALYSISMODULE_H
45 #define GMX_TRAJECTORYANALYSIS_ANALYSISMODULE_H
46
47 #include <string>
48 #include <vector>
49
50 #include "../selection/selection.h" // For gmx::SelectionList
51 #include "../utility/common.h"
52 #include "../utility/uniqueptr.h"
53
54 struct t_pbc;
55 struct t_trxframe;
56
57 namespace gmx
58 {
59
60 class AbstractAnalysisData;
61 class AnalysisData;
62 class AnalysisDataHandle;
63 class AnalysisDataParallelOptions;
64 class Options;
65 class SelectionCollection;
66 class TopologyInformation;
67 class TrajectoryAnalysisModule;
68 class TrajectoryAnalysisSettings;
69
70 /*! \brief
71  * Base class for thread-local data storage during trajectory analysis.
72  *
73  * Thread-local storage of data handles and selections is implemented in this
74  * class; TrajectoryAnalysisModule instances can access the thread-local values
75  * in their TrajectoryAnalysisModule::analyzeFrame() method using dataHandle()
76  * and parallelSelection().
77  *
78  * \see TrajectoryAnalysisModule::startFrames()
79  * \see TrajectoryAnalysisModule::analyzeFrame()
80  * \see TrajectoryAnalysisModule::finishFrames()
81  *
82  * \inpublicapi
83  * \ingroup module_trajectoryanalysis
84  */
85 class TrajectoryAnalysisModuleData
86 {
87     public:
88         virtual ~TrajectoryAnalysisModuleData();
89
90         /*! \brief
91          * Performs any finishing actions after all frames have been processed.
92          *
93          * \throws  unspecified Implementation may throw exceptions to indicate
94          *      errors.
95          *
96          * This function is called immediately before the destructor, after
97          * TrajectoryAnalysisModule::finishFrames().
98          * Derived classes should implement any final operations that need to
99          * be done after successful analysis.
100          * All implementations should call finishDataHandles().
101          */
102         virtual void finish() = 0;
103
104         /*! \brief
105          * Returns a data handle for a given dataset.
106          *
107          * \param[in] data  Analysis data object.
108          * \returns   Data handle for \p data stored in this thread-local data.
109          *
110          * \p data should have previously been registered with
111          * TrajectoryAnalysisModule::registerAnalysisDataset().
112          * If \p data has zero columns in all data sets, the returned data
113          * handle is invalid.
114          *
115          * Does not throw.
116          */
117         AnalysisDataHandle dataHandle(const AnalysisData &data);
118         /*! \brief
119          * Returns a selection that corresponds to the given selection.
120          *
121          * \param[in] selection Global selection object.
122          * \returns   Selection object corresponding to this thread-local data.
123          *
124          * \p selection is the selection object that was obtained from
125          * SelectionOption.  The return value is the corresponding selection
126          * in the selection collection with which this data object was
127          * constructed with.
128          *
129          * Does not throw.
130          */
131         Selection parallelSelection(const Selection &selection);
132         /*! \brief
133          * Returns a set of selection that corresponds to the given selections.
134          *
135          * \throws std::bad_alloc if out of memory.
136          *
137          * Works as parallelSelection(), but for a list of selections at once.
138          *
139          * \see parallelSelection()
140          */
141         SelectionList parallelSelections(const SelectionList &selections);
142
143     protected:
144         /*! \brief
145          * Initializes thread-local storage for data handles and selections.
146          *
147          * \param[in] module     Analysis module to use for data objects.
148          * \param[in] opt        Data parallelization options.
149          * \param[in] selections Thread-local selection collection.
150          * \throws  std::bad_alloc if out of memory.
151          * \throws  unspecified Can throw any exception thrown by
152          *      AnalysisData::startData().
153          *
154          * Calls AnalysisData::startData() on all data objects registered with
155          * TrajectoryAnalysisModule::registerAnalysisDataset() in \p module.
156          * The handles are accessible through dataHandle().
157          */
158         TrajectoryAnalysisModuleData(TrajectoryAnalysisModule          *module,
159                                      const AnalysisDataParallelOptions &opt,
160                                      const SelectionCollection         &selections);
161
162         /*! \brief
163          * Calls finishData() on all data handles.
164          *
165          * \throws  unspecified Can throw any exception thrown by
166          *      AnalysisDataHandle::finishData().
167          *
168          * This function should be called from the implementation of finish()
169          * in all subclasses.
170          */
171         void finishDataHandles();
172
173     private:
174         class Impl;
175
176         PrivateImplPointer<Impl> impl_;
177 };
178
179 //! Smart pointer to manage a TrajectoryAnalysisModuleData object.
180 typedef gmx_unique_ptr<TrajectoryAnalysisModuleData>::type
181     TrajectoryAnalysisModuleDataPointer;
182
183 /*! \brief
184  * Base class for trajectory analysis modules.
185  *
186  * Trajectory analysis methods should derive from this class and override the
187  * necessary virtual methods to implement initialization (initOptions(),
188  * optionsFinished(), initAnalysis(), initAfterFirstFrame()), per-frame analysis
189  * (analyzeFrame()), and final processing (finishFrames(), finishAnalysis(),
190  * writeOutput()).
191  *
192  * For parallel analysis using threads, only a single object is constructed,
193  * but the methods startFrames(), analyzeFrame() and finishFrames() are called
194  * in each thread.  Frame-local data should be initialized in startFrames() and
195  * stored in a class derived from TrajectoryAnalysisModuleData that is passed
196  * to the other methods.  The default implementation of startFrames() can be
197  * used if only data handles and selections need to be thread-local.
198  *
199  * To get the full benefit from this class,
200  * \ref module_analysisdata "analysis data objects" and
201  * \ref module_selection "selections" should be used in the implementation.
202  * See the corresponding modules' documentation for details of how they work.
203  *
204  * Typical way of using AnalysisData in derived classes is to have the
205  * AnalysisData object as a member variable and register it using
206  * registerAnalysisDataset().  Analysis modules are initialized in
207  * initAnalysis() and the processing chain is initialized.  If any of the
208  * modules is required, e.g., for post-processing in finishAnalysis(), it can
209  * be stored in a member variable.  To add data to the data object in
210  * analyzeFrame(), a data handle is obtained using
211  * TrajectoryAnalysisModuleData::dataHandle().
212  *
213  * Typical way of using selections in derived classes is to have the required
214  * \ref Selection objects (or ::SelectionList objects) as member variables, and
215  * add the required selection options in initOptions().  These member variables
216  * can be accessed in initAnalysis() to get general information about the
217  * selections.  In analyzeFrame(), these selection objects should not be used
218  * directly, but instead TrajectoryAnalysisModuleData::parallelSelection()
219  * should be used to obtain a selection object that works correctly also for
220  * parallel analysis.
221  *
222  * Derived classes should use exceptions to indicate errors in the virtual
223  * methods.
224  *
225  * \inpublicapi
226  * \ingroup module_trajectoryanalysis
227  */
228 class TrajectoryAnalysisModule
229 {
230     public:
231         virtual ~TrajectoryAnalysisModule();
232
233         /*! \brief
234          * Initializes options understood by the module.
235          *
236          * \param[in,out] options  Options object to add the options to.
237          * \param[in,out] settings Settings to pass to and from the module.
238          *
239          * This method is called first after the constructor, and it should
240          * add options understood by the module to \p options.  Output values
241          * from options (including selections) should be stored in member
242          * variables.
243          *
244          * In addition to initializing the options, this method can also
245          * provide information about the module's requirements using the
246          * \p settings object; see TrajectoryAnalysisSettings for more details.
247          *
248          * If settings depend on the option values provided by the user, see
249          * optionsFinished().
250          */
251         virtual void initOptions(Options                    *options,
252                                  TrajectoryAnalysisSettings *settings) = 0;
253         /*! \brief
254          * Called after all option values have been set.
255          *
256          * \param[in,out] options  Options object in which options are stored.
257          * \param[in,out] settings Settings to pass to and from the module.
258          *
259          * This method is called after option values have been assigned (but
260          * interactive selection input has not yet been performed).
261          *
262          * If the module needs to change settings that affect topology loading
263          * (can be done using the \p settings object) or selection
264          * initialization (can be done using SelectionOptionInfo) based on
265          * option values, this method has to be overridden.
266          *
267          * The default implementation does nothing.
268          */
269         virtual void optionsFinished(Options                    *options,
270                                      TrajectoryAnalysisSettings *settings);
271         /*! \brief
272          * Initializes the analysis.
273          *
274          * \param[in]    settings Settings to pass to and from the module.
275          * \param[in]    top      Topology information.
276          *
277          * When this function is called, selections have been initialized based
278          * on user input, and a topology has been loaded if provided by the
279          * user.  For dynamic selections, the selections have been evaluated to
280          * the largest possible selection, i.e., the selections passed to
281          * analyzeFrame() are always a subset of the selections provided here.
282          */
283         virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
284                                   const TopologyInformation        &top) = 0;
285         /*! \brief
286          * Performs additional initialization after reading the first frame.
287          *
288          * When this function is called, selections are the same as in
289          * initAnalysis(), i.e., they have not been evaluated for the first
290          * frame.
291          *
292          * It is necessary to override this method only if the module needs to
293          * do initialization for which it requires data from the first frame.
294          *
295          * The default implementation does nothing.
296          */
297         virtual void initAfterFirstFrame(const t_trxframe &fr);
298
299         /*! \brief
300          * Starts the analysis of frames.
301          *
302          * \param[in]  opt
303          * \param[in]  selections  Frame-local selection collection object.
304          * \returns    Data structure for thread-local data.
305          *
306          * This function is necessary only for threaded parallelization.
307          * It is called once for each thread and should initialize a class that
308          * contains any required frame-local data in the returned value.
309          * The default implementation creates a basic data structure that holds
310          * thread-local data handles for all data objects registered with
311          * registerAnalysisDataset(), as well as the thread-local selection
312          * collection.  These can be accessed in analyzeFrame() using the
313          * methods in TrajectoryAnalysisModuleData.
314          * If other thread-local data is needed, this function should be
315          * overridden and it should create an instance of a class derived from
316          * TrajectoryAnalysisModuleData.
317          *
318          * \see TrajectoryAnalysisModuleData
319          */
320         virtual TrajectoryAnalysisModuleDataPointer startFrames(
321             const AnalysisDataParallelOptions &opt,
322             const SelectionCollection         &selections);
323         /*! \brief
324          * Analyzes a single frame.
325          *
326          * \param[in]     frnr   Frame number, a zero-based index that
327          *      uniquely identifies the frame.
328          * \param[in]     fr     Current frame.
329          * \param[in]     pbc    Periodic boundary conditions for \p fr.
330          * \param[in,out] pdata  Data structure for frame-local data.
331          *
332          * This method is called once for each frame to be analyzed, and should
333          * analyze the positions provided in the selections.  Data handles and
334          * selections should be obtained from the \p pdata structure.
335          *
336          * For threaded analysis, this method is called asynchronously in
337          * different threads to analyze different frames.  The \p pdata
338          * structure is one of the structures created with startFrames(),
339          * but no assumptions should be made about which of these data
340          * structures is used.  It is guaranteed that two instances of
341          * analyzeFrame() are not running concurrently with the same \p pdata
342          * data structure.
343          * Any access to data structures not stored in \p pdata should be
344          * designed to be thread-safe.
345          */
346         virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
347                                   TrajectoryAnalysisModuleData *pdata) = 0;
348         /*! \brief
349          * Finishes the analysis of frames.
350          *
351          * \param[in]  pdata    Data structure for thread-local data.
352          *
353          * This method is called once for each call of startFrames(), with the
354          * data structure returned by the corresponding startFrames().
355          * The \p pdata object should be destroyed by the caller after this
356          * function has been called.
357          *
358          * You only need to override this method if you need custom
359          * operations to combine data from the frame-local data structures
360          * to get the final result.  In such cases, the data should be
361          * aggregated in this function and stored in a member attribute.
362          *
363          * The default implementation does nothing.
364          *
365          * \see startFrames()
366          */
367         virtual void finishFrames(TrajectoryAnalysisModuleData *pdata);
368
369         /*! \brief
370          * Postprocesses data after frames have been read.
371          *
372          * \param[in]  nframes  Total number of frames processed.
373          *
374          * This function is called after all finishFrames() calls have been
375          * called.
376          * \p nframes will equal the number of calls to analyzeFrame() that
377          * have occurred.
378          */
379         virtual void finishAnalysis(int nframes) = 0;
380         /*! \brief
381          * Writes output into files and/or standard output/error.
382          *
383          * All output from the module, excluding data written out for each
384          * frame during analyzeFrame(), should be confined into this function.
385          * This function is guaranteed to be called only after
386          * finishAnalysis().
387          */
388         virtual void writeOutput() = 0;
389
390         /*! \brief
391          * Returns the name of the analysis module.
392          *
393          * Does not throw.
394          */
395         const char *name() const;
396         /*! \brief
397          * Returns short description for the analysis module.
398          *
399          * Does not throw.
400          */
401         const char *description() const;
402         /*! \brief
403          * Returns the number of datasets provided by the module.
404          *
405          * Does not throw.
406          */
407         int datasetCount() const;
408         /*! \brief
409          * Returns a vector with the names of datasets provided by the module.
410          *
411          * Does not throw.
412          */
413         const std::vector<std::string> &datasetNames() const;
414         /*! \brief
415          * Returns a pointer to the data set \p index.
416          *
417          * \param[in] index  Data set to query for.
418          * \returns   Reference to the requested data set.
419          * \throws    APIError if \p index is not valid.
420          *
421          * \p index should be >= 0 and < datasetCount().
422          *
423          * The return value is not const to allow callers to add modules to the
424          * data sets. However, the AbstractAnalysisData interface does not
425          * provide any means to alter the data, so the module does not need to
426          * care about external modifications.
427          */
428         AbstractAnalysisData &datasetFromIndex(int index) const;
429         /*! \brief
430          * Returns a pointer to the data set with name \p name
431          *
432          * \param[in] name  Data set to query for.
433          * \returns   Reference to the requested data set.
434          * \throws    APIError if \p name is not valid.
435          *
436          * \p name should be one of the names returned by datasetNames().
437          *
438          * The return value is not const to allow callers to add modules to the
439          * data sets. However, the AbstractAnalysisData interface does not
440          * provide any means to alter the data, so the module does not need to
441          * care about external modifications.
442          */
443         AbstractAnalysisData &datasetFromName(const char *name) const;
444
445     protected:
446         /*! \brief
447          * Initializes the dataset registration mechanism.
448          *
449          * \param[in] name         Name for the module.
450          * \param[in] description  One-line description for the module.
451          * \throws    std::bad_alloc if out of memory.
452          */
453         TrajectoryAnalysisModule(const char *name, const char *description);
454
455         /*! \brief
456          * Registers a dataset that exports data.
457          *
458          * \param     data  Data object to register.
459          * \param[in] name  Name to register the dataset with.
460          * \throws    std::bad_alloc if out of memory.
461          *
462          * Registers \p data as a dataset that provides output from the
463          * analysis module.  Callers for the module can access the dataset
464          * with datasetFromName() using \p name as an AbstractAnalysisData
465          * object.  This allows them to add their own data modules to do extra
466          * processing.
467          *
468          * \p name must be unique across all calls within the same
469          * TrajectoryAnalysisModule instance.
470          */
471         void registerBasicDataset(AbstractAnalysisData *data, const char *name);
472         /*! \brief
473          * Registers a parallelized dataset that exports data.
474          *
475          * \param     data  AnalysisData object to register.
476          * \param[in] name  Name to register the dataset with.
477          * \throws    std::bad_alloc if out of memory.
478          *
479          * This method works as registerBasicDataset(), but additionally allows
480          * data handles for \p data to be accessed using
481          * TrajectoryAnalysisData.
482          *
483          * \see registerBasicDataset()
484          */
485         void registerAnalysisDataset(AnalysisData *data, const char *name);
486
487     private:
488         class Impl;
489
490         PrivateImplPointer<Impl> impl_;
491
492         /*! \brief
493          * Needed to access the registered analysis data sets.
494          */
495         friend class TrajectoryAnalysisModuleData;
496 };
497
498 //! Smart pointer to manage a TrajectoryAnalysisModule.
499 typedef gmx_unique_ptr<TrajectoryAnalysisModule>::type
500     TrajectoryAnalysisModulePointer;
501
502 } // namespace gmx
503
504 #endif