Add basic grouping to Options
[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,2015, 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 <boost/shared_ptr.hpp>
51
52 #include "gromacs/selection/selection.h" // For gmx::SelectionList
53 #include "gromacs/utility/classhelpers.h"
54
55 struct t_pbc;
56 struct t_trxframe;
57
58 namespace gmx
59 {
60
61 class AbstractAnalysisData;
62 class AnalysisData;
63 class AnalysisDataHandle;
64 class AnalysisDataParallelOptions;
65 class IOptionsContainer;
66 class Options;
67 class SelectionCollection;
68 class TopologyInformation;
69 class TrajectoryAnalysisModule;
70 class TrajectoryAnalysisSettings;
71
72 /*! \brief
73  * Base class for thread-local data storage during trajectory analysis.
74  *
75  * Thread-local storage of data handles and selections is implemented in this
76  * class; TrajectoryAnalysisModule instances can access the thread-local values
77  * in their TrajectoryAnalysisModule::analyzeFrame() method using dataHandle()
78  * and parallelSelection().
79  *
80  * \see TrajectoryAnalysisModule::startFrames()
81  * \see TrajectoryAnalysisModule::analyzeFrame()
82  * \see TrajectoryAnalysisModule::finishFrames()
83  *
84  * \inpublicapi
85  * \ingroup module_trajectoryanalysis
86  */
87 class TrajectoryAnalysisModuleData
88 {
89     public:
90         virtual ~TrajectoryAnalysisModuleData();
91
92         /*! \brief
93          * Performs any finishing actions after all frames have been processed.
94          *
95          * \throws  unspecified Implementation may throw exceptions to indicate
96          *      errors.
97          *
98          * This function is called immediately before the destructor, after
99          * TrajectoryAnalysisModule::finishFrames().
100          * Derived classes should implement any final operations that need to
101          * be done after successful analysis.
102          * All implementations should call finishDataHandles().
103          */
104         virtual void finish() = 0;
105
106         /*! \brief
107          * Returns a data handle for a given dataset.
108          *
109          * \param[in] data  Analysis data object.
110          * \returns   Data handle for \p data stored in this thread-local data.
111          *
112          * \p data should have previously been registered with
113          * TrajectoryAnalysisModule::registerAnalysisDataset().
114          * If \p data has zero columns in all data sets, the returned data
115          * handle is invalid.
116          *
117          * Does not throw.
118          */
119         AnalysisDataHandle dataHandle(const AnalysisData &data);
120         /*! \brief
121          * Returns a selection that corresponds to the given selection.
122          *
123          * \param[in] selection Global selection object.
124          * \returns   Selection object corresponding to this thread-local data.
125          *
126          * \p selection is the selection object that was obtained from
127          * SelectionOption.  The return value is the corresponding selection
128          * in the selection collection with which this data object was
129          * constructed with.
130          *
131          * Does not throw.
132          */
133         Selection parallelSelection(const Selection &selection);
134         /*! \brief
135          * Returns a set of selection that corresponds to the given selections.
136          *
137          * \throws std::bad_alloc if out of memory.
138          *
139          * Works as parallelSelection(), but for a list of selections at once.
140          *
141          * \see parallelSelection()
142          */
143         SelectionList parallelSelections(const SelectionList &selections);
144
145     protected:
146         /*! \brief
147          * Initializes thread-local storage for data handles and selections.
148          *
149          * \param[in] module     Analysis module to use for data objects.
150          * \param[in] opt        Data parallelization options.
151          * \param[in] selections Thread-local selection collection.
152          * \throws  std::bad_alloc if out of memory.
153          * \throws  unspecified Can throw any exception thrown by
154          *      AnalysisData::startData().
155          *
156          * Calls AnalysisData::startData() on all data objects registered with
157          * TrajectoryAnalysisModule::registerAnalysisDataset() in \p module.
158          * The handles are accessible through dataHandle().
159          */
160         TrajectoryAnalysisModuleData(TrajectoryAnalysisModule          *module,
161                                      const AnalysisDataParallelOptions &opt,
162                                      const SelectionCollection         &selections);
163
164         /*! \brief
165          * Calls finishData() on all data handles.
166          *
167          * \throws  unspecified Can throw any exception thrown by
168          *      AnalysisDataHandle::finishData().
169          *
170          * This function should be called from the implementation of finish()
171          * in all subclasses.
172          */
173         void finishDataHandles();
174
175     private:
176         class Impl;
177
178         PrivateImplPointer<Impl> impl_;
179 };
180
181 //! Smart pointer to manage a TrajectoryAnalysisModuleData object.
182 typedef boost::shared_ptr<TrajectoryAnalysisModuleData>
183     TrajectoryAnalysisModuleDataPointer;
184
185 /*! \brief
186  * Base class for trajectory analysis modules.
187  *
188  * Trajectory analysis methods should derive from this class and override the
189  * necessary virtual methods to implement initialization (initOptions(),
190  * optionsFinished(), initAnalysis(), initAfterFirstFrame()), per-frame analysis
191  * (analyzeFrame()), and final processing (finishFrames(), finishAnalysis(),
192  * writeOutput()).
193  *
194  * For parallel analysis using threads, only a single object is constructed,
195  * but the methods startFrames(), analyzeFrame() and finishFrames() are called
196  * in each thread.  Frame-local data should be initialized in startFrames() and
197  * stored in a class derived from TrajectoryAnalysisModuleData that is passed
198  * to the other methods.  The default implementation of startFrames() can be
199  * used if only data handles and selections need to be thread-local.
200  *
201  * To get the full benefit from this class,
202  * \ref module_analysisdata "analysis data objects" and
203  * \ref module_selection "selections" should be used in the implementation.
204  * See the corresponding modules' documentation for details of how they work.
205  *
206  * Typical way of using AnalysisData in derived classes is to have the
207  * AnalysisData object as a member variable and register it using
208  * registerAnalysisDataset().  Analysis modules are initialized in
209  * initAnalysis() and the processing chain is initialized.  If any of the
210  * modules is required, e.g., for post-processing in finishAnalysis(), it can
211  * be stored in a member variable.  To add data to the data object in
212  * analyzeFrame(), a data handle is obtained using
213  * TrajectoryAnalysisModuleData::dataHandle().
214  *
215  * Typical way of using selections in derived classes is to have the required
216  * \ref Selection objects (or ::SelectionList objects) as member variables, and
217  * add the required selection options in initOptions().  These member variables
218  * can be accessed in initAnalysis() to get general information about the
219  * selections.  In analyzeFrame(), these selection objects should not be used
220  * directly, but instead TrajectoryAnalysisModuleData::parallelSelection()
221  * should be used to obtain a selection object that works correctly also for
222  * parallel analysis.
223  *
224  * Derived classes should use exceptions to indicate errors in the virtual
225  * methods.
226  *
227  * \inpublicapi
228  * \ingroup module_trajectoryanalysis
229  */
230 class TrajectoryAnalysisModule
231 {
232     public:
233         virtual ~TrajectoryAnalysisModule();
234
235         /*! \brief
236          * Initializes options understood by the module.
237          *
238          * \param[in,out] options  Options object to add the options to.
239          * \param[in,out] settings Settings to pass to and from the module.
240          *
241          * This method is called first after the constructor, and it should
242          * add options understood by the module to \p options.  Output values
243          * from options (including selections) should be stored in member
244          * variables.
245          *
246          * In addition to initializing the options, this method can also
247          * provide information about the module's requirements using the
248          * \p settings object; see TrajectoryAnalysisSettings for more details.
249          *
250          * If settings depend on the option values provided by the user, see
251          * optionsFinished().
252          */
253         virtual void initOptions(IOptionsContainer          *options,
254                                  TrajectoryAnalysisSettings *settings) = 0;
255         /*! \brief
256          * Called after all option values have been set.
257          *
258          * \param[in,out] options  Options object in which options are stored.
259          * \param[in,out] settings Settings to pass to and from the module.
260          *
261          * This method is called after option values have been assigned (but
262          * interactive selection input has not yet been performed).
263          *
264          * If the module needs to change settings that affect topology loading
265          * (can be done using the \p settings object) or selection
266          * initialization (can be done using SelectionOptionInfo) based on
267          * option values, this method has to be overridden.
268          *
269          * The default implementation does nothing.
270          */
271         virtual void optionsFinished(Options                    *options,
272                                      TrajectoryAnalysisSettings *settings);
273         /*! \brief
274          * Initializes the analysis.
275          *
276          * \param[in]    settings Settings to pass to and from the module.
277          * \param[in]    top      Topology information.
278          *
279          * When this function is called, selections have been initialized based
280          * on user input, and a topology has been loaded if provided by the
281          * user.  For dynamic selections, the selections have been evaluated to
282          * the largest possible selection, i.e., the selections passed to
283          * analyzeFrame() are always a subset of the selections provided here.
284          */
285         virtual void initAnalysis(const TrajectoryAnalysisSettings &settings,
286                                   const TopologyInformation        &top) = 0;
287         /*! \brief
288          * Performs additional initialization after reading the first frame.
289          *
290          * When this function is called, selections are the same as in
291          * initAnalysis(), i.e., they have not been evaluated for the first
292          * frame.
293          *
294          * It is necessary to override this method only if the module needs to
295          * do initialization for which it requires data from the first frame.
296          *
297          * The default implementation does nothing.
298          */
299         virtual void initAfterFirstFrame(const TrajectoryAnalysisSettings &settings,
300                                          const t_trxframe                 &fr);
301
302         /*! \brief
303          * Starts the analysis of frames.
304          *
305          * \param[in]  opt
306          * \param[in]  selections  Frame-local selection collection object.
307          * \returns    Data structure for thread-local data.
308          *
309          * This function is necessary only for threaded parallelization.
310          * It is called once for each thread and should initialize a class that
311          * contains any required frame-local data in the returned value.
312          * The default implementation creates a basic data structure that holds
313          * thread-local data handles for all data objects registered with
314          * registerAnalysisDataset(), as well as the thread-local selection
315          * collection.  These can be accessed in analyzeFrame() using the
316          * methods in TrajectoryAnalysisModuleData.
317          * If other thread-local data is needed, this function should be
318          * overridden and it should create an instance of a class derived from
319          * TrajectoryAnalysisModuleData.
320          *
321          * \see TrajectoryAnalysisModuleData
322          */
323         virtual TrajectoryAnalysisModuleDataPointer startFrames(
324             const AnalysisDataParallelOptions &opt,
325             const SelectionCollection         &selections);
326         /*! \brief
327          * Analyzes a single frame.
328          *
329          * \param[in]     frnr   Frame number, a zero-based index that
330          *      uniquely identifies the frame.
331          * \param[in]     fr     Current frame.
332          * \param[in]     pbc    Periodic boundary conditions for \p fr.
333          * \param[in,out] pdata  Data structure for frame-local data.
334          *
335          * This method is called once for each frame to be analyzed, and should
336          * analyze the positions provided in the selections.  Data handles and
337          * selections should be obtained from the \p pdata structure.
338          *
339          * For threaded analysis, this method is called asynchronously in
340          * different threads to analyze different frames.  The \p pdata
341          * structure is one of the structures created with startFrames(),
342          * but no assumptions should be made about which of these data
343          * structures is used.  It is guaranteed that two instances of
344          * analyzeFrame() are not running concurrently with the same \p pdata
345          * data structure.
346          * Any access to data structures not stored in \p pdata should be
347          * designed to be thread-safe.
348          */
349         virtual void analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
350                                   TrajectoryAnalysisModuleData *pdata) = 0;
351         /*! \brief
352          * Finishes the analysis of frames.
353          *
354          * \param[in]  pdata    Data structure for thread-local data.
355          *
356          * This method is called once for each call of startFrames(), with the
357          * data structure returned by the corresponding startFrames().
358          * The \p pdata object should be destroyed by the caller after this
359          * function has been called.
360          *
361          * You only need to override this method if you need custom
362          * operations to combine data from the frame-local data structures
363          * to get the final result.  In such cases, the data should be
364          * aggregated in this function and stored in a member attribute.
365          *
366          * The default implementation does nothing.
367          *
368          * \see startFrames()
369          */
370         virtual void finishFrames(TrajectoryAnalysisModuleData *pdata);
371
372         /*! \brief
373          * Postprocesses data after frames have been read.
374          *
375          * \param[in]  nframes  Total number of frames processed.
376          *
377          * This function is called after all finishFrames() calls have been
378          * called.
379          * \p nframes will equal the number of calls to analyzeFrame() that
380          * have occurred.
381          */
382         virtual void finishAnalysis(int nframes) = 0;
383         /*! \brief
384          * Writes output into files and/or standard output/error.
385          *
386          * All output from the module, excluding data written out for each
387          * frame during analyzeFrame(), should be confined into this function.
388          * This function is guaranteed to be called only after
389          * finishAnalysis().
390          */
391         virtual void writeOutput() = 0;
392
393         /*! \brief
394          * Returns the name of the analysis module.
395          *
396          * Does not throw.
397          */
398         const char *name() const;
399         /*! \brief
400          * Returns short description for the analysis module.
401          *
402          * Does not throw.
403          */
404         const char *description() const;
405         /*! \brief
406          * Returns the number of datasets provided by the module.
407          *
408          * Does not throw.
409          */
410         int datasetCount() const;
411         /*! \brief
412          * Returns a vector with the names of datasets provided by the module.
413          *
414          * Does not throw.
415          */
416         const std::vector<std::string> &datasetNames() const;
417         /*! \brief
418          * Returns a pointer to the data set \p index.
419          *
420          * \param[in] index  Data set to query for.
421          * \returns   Reference to the requested data set.
422          * \throws    APIError if \p index is not valid.
423          *
424          * \p index should be >= 0 and < datasetCount().
425          *
426          * The return value is not const to allow callers to add modules to the
427          * data sets. However, the AbstractAnalysisData interface does not
428          * provide any means to alter the data, so the module does not need to
429          * care about external modifications.
430          */
431         AbstractAnalysisData &datasetFromIndex(int index) const;
432         /*! \brief
433          * Returns a pointer to the data set with name \p name
434          *
435          * \param[in] name  Data set to query for.
436          * \returns   Reference to the requested data set.
437          * \throws    APIError if \p name is not valid.
438          *
439          * \p name should be one of the names returned by datasetNames().
440          *
441          * The return value is not const to allow callers to add modules to the
442          * data sets. However, the AbstractAnalysisData interface does not
443          * provide any means to alter the data, so the module does not need to
444          * care about external modifications.
445          */
446         AbstractAnalysisData &datasetFromName(const char *name) const;
447         /*! \brief
448          * Processes data in AnalysisData objects in serial for each frame.
449          *
450          * \param[in] frameIndex  Index of the frame that has been finished.
451          *
452          * This method is called by the framework in order for each frame,
453          * after the analysis for that frame has been finished.  These calls
454          * always execute in serial and in sequential frame order, even during
455          * parallel analysis where multiple analyzeFrame() calls may be
456          * executing concurrently.
457          *
458          * \see AnalysisData::finishFrameSerial()
459          */
460         void finishFrameSerial(int frameIndex);
461
462     protected:
463         /*! \brief
464          * Initializes the dataset registration mechanism.
465          *
466          * \param[in] name         Name for the module.
467          * \param[in] description  One-line description for the module.
468          * \throws    std::bad_alloc if out of memory.
469          */
470         TrajectoryAnalysisModule(const char *name, const char *description);
471
472         /*! \brief
473          * Registers a dataset that exports data.
474          *
475          * \param     data  Data object to register.
476          * \param[in] name  Name to register the dataset with.
477          * \throws    std::bad_alloc if out of memory.
478          *
479          * Registers \p data as a dataset that provides output from the
480          * analysis module.  Callers for the module can access the dataset
481          * with datasetFromName() using \p name as an AbstractAnalysisData
482          * object.  This allows them to add their own data modules to do extra
483          * processing.
484          *
485          * \p name must be unique across all calls within the same
486          * TrajectoryAnalysisModule instance.
487          */
488         void registerBasicDataset(AbstractAnalysisData *data, const char *name);
489         /*! \brief
490          * Registers a parallelized dataset that exports data.
491          *
492          * \param     data  AnalysisData object to register.
493          * \param[in] name  Name to register the dataset with.
494          * \throws    std::bad_alloc if out of memory.
495          *
496          * This method works as registerBasicDataset(), but additionally allows
497          * data handles for \p data to be accessed using
498          * TrajectoryAnalysisData.
499          *
500          * \see registerBasicDataset()
501          */
502         void registerAnalysisDataset(AnalysisData *data, const char *name);
503
504     private:
505         class Impl;
506
507         PrivateImplPointer<Impl> impl_;
508
509         /*! \brief
510          * Needed to access the registered analysis data sets.
511          */
512         friend class TrajectoryAnalysisModuleData;
513 };
514
515 //! Smart pointer to manage a TrajectoryAnalysisModule.
516 typedef boost::shared_ptr<TrajectoryAnalysisModule>
517     TrajectoryAnalysisModulePointer;
518
519 } // namespace gmx
520
521 #endif