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