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