7462f19badb02604ef287afe0524015bc4c2d57e
[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.
5  * Copyright (c) 2015,2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \file
37  * \brief
38  * Declares gmx::TrajectoryAnalysisModule and
39  * gmx::TrajectoryAnalysisModuleData.
40  *
41  * \author Teemu Murtola <teemu.murtola@gmail.com>
42  * \inpublicapi
43  * \ingroup module_trajectoryanalysis
44  */
45 #ifndef GMX_TRAJECTORYANALYSIS_ANALYSISMODULE_H
46 #define GMX_TRAJECTORYANALYSIS_ANALYSISMODULE_H
47
48 #include <memory>
49 #include <string>
50 #include <vector>
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     static 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     static 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 std::unique_ptr<TrajectoryAnalysisModuleData> 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(IOptionsContainer* options, TrajectoryAnalysisSettings* settings) = 0;
253     /*! \brief
254      * Called after all option values have been set.
255      *
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(TrajectoryAnalysisSettings* settings);
269     /*! \brief
270      * Initializes the analysis.
271      *
272      * \param[in]    settings Settings to pass to and from the module.
273      * \param[in]    top      Topology information.
274      *
275      * When this function is called, selections have been initialized based
276      * on user input, and a topology has been loaded if provided by the
277      * user.  For dynamic selections, the selections have been evaluated to
278      * the largest possible selection, i.e., the selections passed to
279      * analyzeFrame() are always a subset of the selections provided here.
280      */
281     virtual void initAnalysis(const TrajectoryAnalysisSettings& settings,
282                               const TopologyInformation&        top) = 0;
283     /*! \brief
284      * Performs additional initialization after reading the first frame.
285      *
286      * When this function is called, selections are the same as in
287      * initAnalysis(), i.e., they have not been evaluated for the first
288      * frame.
289      *
290      * It is necessary to override this method only if the module needs to
291      * do initialization for which it requires data from the first frame.
292      *
293      * The default implementation does nothing.
294      */
295     virtual void initAfterFirstFrame(const TrajectoryAnalysisSettings& settings, const t_trxframe& fr);
296
297     /*! \brief
298      * Starts the analysis of frames.
299      *
300      * \param[in]  opt         Parallel options
301      * \param[in]  selections  Frame-local selection collection object.
302      * \returns    Data structure for thread-local data.
303      *
304      * This function is necessary only for threaded parallelization.
305      * It is called once for each thread and should initialize a class that
306      * contains any required frame-local data in the returned value.
307      * The default implementation creates a basic data structure that holds
308      * thread-local data handles for all data objects registered with
309      * registerAnalysisDataset(), as well as the thread-local selection
310      * collection.  These can be accessed in analyzeFrame() using the
311      * methods in TrajectoryAnalysisModuleData.
312      * If other thread-local data is needed, this function should be
313      * overridden and it should create an instance of a class derived from
314      * TrajectoryAnalysisModuleData.
315      *
316      * \see TrajectoryAnalysisModuleData
317      */
318     virtual TrajectoryAnalysisModuleDataPointer startFrames(const AnalysisDataParallelOptions& opt,
319                                                             const SelectionCollection& selections);
320     /*! \brief
321      * Analyzes a single frame.
322      *
323      * \param[in]     frnr   Frame number, a zero-based index that
324      *      uniquely identifies the frame.
325      * \param[in]     fr     Current frame.
326      * \param[in]     pbc    Periodic boundary conditions for \p fr.
327      * \param[in,out] pdata  Data structure for frame-local data.
328      *
329      * This method is called once for each frame to be analyzed, and should
330      * analyze the positions provided in the selections.  Data handles and
331      * selections should be obtained from the \p pdata structure.
332      *
333      * For threaded analysis, this method is called asynchronously in
334      * different threads to analyze different frames.  The \p pdata
335      * structure is one of the structures created with startFrames(),
336      * but no assumptions should be made about which of these data
337      * structures is used.  It is guaranteed that two instances of
338      * analyzeFrame() are not running concurrently with the same \p pdata
339      * data structure.
340      * Any access to data structures not stored in \p pdata should be
341      * designed to be thread-safe.
342      */
343     virtual void analyzeFrame(int frnr, const t_trxframe& fr, t_pbc* pbc, TrajectoryAnalysisModuleData* pdata) = 0;
344     /*! \brief
345      * Finishes the analysis of frames.
346      *
347      * \param[in]  pdata    Data structure for thread-local data.
348      *
349      * This method is called once for each call of startFrames(), with the
350      * data structure returned by the corresponding startFrames().
351      * The \p pdata object should be destroyed by the caller after this
352      * function has been called.
353      *
354      * You only need to override this method if you need custom
355      * operations to combine data from the frame-local data structures
356      * to get the final result.  In such cases, the data should be
357      * aggregated in this function and stored in a member attribute.
358      *
359      * The default implementation does nothing.
360      *
361      * \see startFrames()
362      */
363     virtual void finishFrames(TrajectoryAnalysisModuleData* pdata);
364
365     /*! \brief
366      * Postprocesses data after frames have been read.
367      *
368      * \param[in]  nframes  Total number of frames processed.
369      *
370      * This function is called after all finishFrames() calls have been
371      * called.
372      * \p nframes will equal the number of calls to analyzeFrame() that
373      * have occurred.
374      */
375     virtual void finishAnalysis(int nframes) = 0;
376     /*! \brief
377      * Writes output into files and/or standard output/error.
378      *
379      * All output from the module, excluding data written out for each
380      * frame during analyzeFrame(), should be confined into this function.
381      * This function is guaranteed to be called only after
382      * finishAnalysis().
383      */
384     virtual void writeOutput() = 0;
385
386     /*! \brief
387      * Returns the number of datasets provided by the module.
388      *
389      * Does not throw.
390      */
391     int datasetCount() const;
392     /*! \brief
393      * Returns a vector with the names of datasets provided by the module.
394      *
395      * Does not throw.
396      */
397     const std::vector<std::string>& datasetNames() const;
398     /*! \brief
399      * Returns a pointer to the data set \p index.
400      *
401      * \param[in] index  Data set to query for.
402      * \returns   Reference to the requested data set.
403      * \throws    APIError if \p index is not valid.
404      *
405      * \p index should be >= 0 and < datasetCount().
406      *
407      * The return value is not const to allow callers to add modules to the
408      * data sets. However, the AbstractAnalysisData interface does not
409      * provide any means to alter the data, so the module does not need to
410      * care about external modifications.
411      */
412     AbstractAnalysisData& datasetFromIndex(int index) const;
413     /*! \brief
414      * Returns a pointer to the data set with name \p name
415      *
416      * \param[in] name  Data set to query for.
417      * \returns   Reference to the requested data set.
418      * \throws    APIError if \p name is not valid.
419      *
420      * \p name should be one of the names returned by datasetNames().
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& datasetFromName(const char* name) const;
428     /*! \brief
429      * Processes data in AnalysisData objects in serial for each frame.
430      *
431      * \param[in] frameIndex  Index of the frame that has been finished.
432      *
433      * This method is called by the framework in order for each frame,
434      * after the analysis for that frame has been finished.  These calls
435      * always execute in serial and in sequential frame order, even during
436      * parallel analysis where multiple analyzeFrame() calls may be
437      * executing concurrently.
438      *
439      * \see AnalysisData::finishFrameSerial()
440      */
441     void finishFrameSerial(int frameIndex);
442
443 protected:
444     /*! \brief
445      * Initializes the dataset registration mechanism.
446      *
447      * \throws    std::bad_alloc if out of memory.
448      */
449     TrajectoryAnalysisModule();
450
451     /*! \brief
452      * Registers a dataset that exports data.
453      *
454      * \param     data  Data object to register.
455      * \param[in] name  Name to register the dataset with.
456      * \throws    std::bad_alloc if out of memory.
457      *
458      * Registers \p data as a dataset that provides output from the
459      * analysis module.  Callers for the module can access the dataset
460      * with datasetFromName() using \p name as an AbstractAnalysisData
461      * object.  This allows them to add their own data modules to do extra
462      * processing.
463      *
464      * \p name must be unique across all calls within the same
465      * TrajectoryAnalysisModule instance.
466      */
467     void registerBasicDataset(AbstractAnalysisData* data, const char* name);
468     /*! \brief
469      * Registers a parallelized dataset that exports data.
470      *
471      * \param     data  AnalysisData object to register.
472      * \param[in] name  Name to register the dataset with.
473      * \throws    std::bad_alloc if out of memory.
474      *
475      * This method works as registerBasicDataset(), but additionally allows
476      * data handles for \p data to be accessed using
477      * TrajectoryAnalysisData.
478      *
479      * \see registerBasicDataset()
480      */
481     void registerAnalysisDataset(AnalysisData* data, const char* name);
482
483 private:
484     class Impl;
485
486     PrivateImplPointer<Impl> impl_;
487
488     /*! \brief
489      * Needed to access the registered analysis data sets.
490      */
491     friend class TrajectoryAnalysisModuleData;
492 };
493
494 //! Smart pointer to manage a TrajectoryAnalysisModule.
495 typedef std::unique_ptr<TrajectoryAnalysisModule> TrajectoryAnalysisModulePointer;
496
497 } // namespace gmx
498
499 #endif