Merge branch release-4-6
[alexxy/gromacs.git] / src / gromacs / analysisdata / datastorage.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013, by the GROMACS development team, led by
5  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
6  * others, as listed in the AUTHORS file in the top-level source
7  * 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 /*! \libinternal \file
36  * \brief
37  * Declares gmx::AnalysisDataStorage.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \inlibraryapi
41  * \ingroup module_analysisdata
42  */
43 #ifndef GMX_ANALYSISDATA_DATASTORAGE_H
44 #define GMX_ANALYSISDATA_DATASTORAGE_H
45
46 #include <vector>
47
48 #include "../legacyheaders/types/simple.h"
49
50 #include "../utility/common.h"
51 #include "../utility/gmxassert.h"
52
53 #include "dataframe.h"
54
55 namespace gmx
56 {
57
58 class AbstractAnalysisData;
59 class AnalysisDataFrameHeader;
60 class AnalysisDataFrameRef;
61 class AnalysisDataModuleManager;
62 class AnalysisDataParallelOptions;
63
64 class AnalysisDataStorage;
65
66 namespace internal
67 {
68 class AnalysisDataStorageImpl;
69 class AnalysisDataStorageFrameData;
70 }   // namespace internal
71
72 /*! \libinternal \brief
73  * Allows assigning values for a data frame in AnalysisDataStorage.
74  *
75  * This class implements the necessary methods to add new data into the
76  * storage.  AnalysisDataStorage::startFrame() returns an object of this type,
77  * which can be used to add one or more point sets to that data frame.
78  * When all data has been added, finishFrame() needs to be called.
79  *
80  * \inlibraryapi
81  * \ingroup module_analysisdata
82  */
83 class AnalysisDataStorageFrame
84 {
85     public:
86         /*! \brief Frees the frame object.
87          *
88          * Should not be called outside AnalysisDataStorage.
89          */
90         ~AnalysisDataStorageFrame();
91
92         /*! \brief
93          * Select data set that all other methods operate on.
94          *
95          * \param[in] index  Zero-based data set index to select.
96          *
97          * With multipoint data, a single point set can only contain values in
98          * a single data set.
99          * With non-multipoint data, arbitrary sequences of selectDataSet() and
100          * setValue() are supported.  The full frame is notified to the modules
101          * once it is finished.
102          *
103          * Does not throw.
104          */
105         void selectDataSet(int index);
106
107         //! Returns number of columns for the frame.
108         int columnCount() const { return columnCount_; }
109
110         /*! \brief
111          * Sets value for a column.
112          *
113          * \param[in] column  Zero-based column index.
114          * \param[in] value   Value to set for the column.
115          * \param[in] bPresent Present flag to set for the column.
116          *
117          * If called multiple times for a column (within one point set for
118          * multipoint data), old values are overwritten.
119          *
120          * Does not throw.
121          */
122         void setValue(int column, real value, bool bPresent = true)
123         {
124             GMX_ASSERT(column >= 0 && column < columnCount(),
125                        "Invalid column index");
126             values_[currentOffset_ + column].setValue(value, bPresent);
127             bPointSetInProgress_ = true;
128         }
129         /*! \brief
130          * Sets value for a column.
131          *
132          * \param[in] column  Zero-based column index.
133          * \param[in] value   Value to set for the column.
134          * \param[in] error   Error estimate to set for the column.
135          * \param[in] bPresent Present flag to set for the column.
136          *
137          * If called multiple times for a column (within one point set for
138          * multipoint data), old values are overwritten.
139          *
140          * Does not throw.
141          */
142         void setValue(int column, real value, real error, bool bPresent = true)
143         {
144             GMX_ASSERT(column >= 0 && column < columnCount(),
145                        "Invalid column index");
146             values_[currentOffset_ + column].setValue(value, error, bPresent);
147             bPointSetInProgress_ = true;
148         }
149         /*! \brief
150          * Access value for a column.
151          *
152          * \param[in] column  Zero-based column index.
153          *
154          * Should only be called after the column value has been set using
155          * setValue(); assigning a value to \c value(i) does not mark the
156          * column as set.
157          *
158          * Does not throw.
159          */
160         real &value(int column)
161         {
162             GMX_ASSERT(column >= 0 && column < columnCount(),
163                        "Invalid column index");
164             return values_[currentOffset_ + column].value();
165         }
166         /*! \brief
167          * Access value for a column.
168          *
169          * \param[in] column  Zero-based column index.
170          *
171          * Should only be called after the column value has been set using
172          * setValue().
173          *
174          * Does not throw.
175          */
176         real value(int column) const
177         {
178             GMX_ASSERT(column >= 0 && column < columnCount(),
179                        "Invalid column index");
180             return values_[currentOffset_ + column].value();
181         }
182         /*! \brief
183          * Mark point set as finished for multipoint data.
184          *
185          * Must be called after each point set for multipoint data, including
186          * the last (i.e., no values must be set between the last call to this
187          * method and AnalysisDataStorage::finishFrame()).
188          * Must not be called for non-multipoint data.
189          *
190          * After this method has been called, all values appear as not set.
191          *
192          * May call AnalysisDataModuleManager::notifyPointsAdd() and
193          * AnalysisDataModuleManager::notifyParallelPointsAdd(), and may throw
194          * any exception these methods throw.
195          */
196         void finishPointSet();
197         /*! \brief
198          * Finish storing a frame.
199          *
200          * Must be called exactly once for each frame returned by startFrame(),
201          * after the corresponding call.
202          * The frame object must not be accessed after the call.
203          *
204          * Calls notification methods in AnalysisDataModuleManager, and may
205          * throw any exceptions these methods throw.
206          */
207         void finishFrame();
208
209     private:
210
211         /*! \brief
212          * Create a new storage frame.
213          *
214          * \param[in] data  Data object for which the frame is for
215          *      (used for data set and column counts).
216          */
217         explicit AnalysisDataStorageFrame(const AbstractAnalysisData &data);
218
219         //! Clear all column values from the frame.
220         void clearValues();
221
222         //! Implementation data.
223         internal::AnalysisDataStorageFrameData *data_;
224         //! Values for the currently in-progress point set.
225         std::vector<AnalysisDataValue>          values_;
226
227         //! Index of the currently active dataset.
228         int                                     currentDataSet_;
229         //! Offset of the first value in \a values_ for the current data set.
230         int                                     currentOffset_;
231         //! Number of columns in the current data set.
232         int                                     columnCount_;
233
234         //! Whether any values have been set in the current point set.
235         bool                                    bPointSetInProgress_;
236
237         //! Needed for access to the constructor.
238         friend class internal::AnalysisDataStorageImpl;
239         //! Needed for managing the frame the object points to.
240         friend class internal::AnalysisDataStorageFrameData;
241
242         GMX_DISALLOW_COPY_AND_ASSIGN(AnalysisDataStorageFrame);
243 };
244
245 /*! \libinternal \brief
246  * Helper class that implements storage of data.
247  *
248  * This class implements a standard way of storing data to avoid implementing
249  * storage in each class derived from AbstractAnalysisData separately.
250  * To use this class in a class derived from AbstractAnalysisData, a member
251  * variable of this type should be declared and the pure virtual methods
252  * forwarded to frameCount(), tryGetDataFrame() and requestStorage().
253  * Storage properties should be set up, and then startDataStorage() or
254  * startParallelDataStorage() called.
255  * New frames can then be added using startFrame(), currentFrame() and
256  * finishFrame() methods.  When all frames are ready, finishDataStorage() must
257  * be called.  These methods (and AnalysisDataStorageFrame::finishPointSet())
258  * take the responsibility of calling all the notification methods in
259  * AnalysisDataModuleManager,
260  *
261  * \todo
262  * Proper multi-threaded implementation.
263  *
264  * \inlibraryapi
265  * \ingroup module_analysisdata
266  */
267 class AnalysisDataStorage
268 {
269     public:
270         //! Constructs a storage object.
271         AnalysisDataStorage();
272         ~AnalysisDataStorage();
273
274         /*! \brief
275          * Returns the number of ready frames.
276          *
277          * This method is designed such that calls to
278          * AbstractAnalysisData::frameCount() can be directly forwarded to this
279          * method.  See that method for more documentation.
280          *
281          * If this method returns N, this means that the first N frames have
282          * all been finished.
283          *
284          * \see AbstractAnalysisData::frameCount()
285          */
286         int frameCount() const;
287         /*! \brief
288          * Implements access to data frames.
289          *
290          * This method is designed such that calls to
291          * AbstractAnalysisData::tryGetDataFrameInternal() can be directly
292          * forwarded to this method.  See that method for more documentation.
293          *
294          * A valid reference for a frame will be returned after finishFrame()
295          * has been called for that frame.
296          *
297          * \see AbstractAnalysisData::tryGetDataFrameInternal()
298          */
299         AnalysisDataFrameRef tryGetDataFrame(int index) const;
300         /*! \brief
301          * Implements storage requests.
302          *
303          * This method is designed such that calls to
304          * AbstractAnalysisData::requestStorageInternal() can be directly
305          * forwarded to this method.  See that method for more documentation.
306          *
307          * \see AbstractAnalysisData::requestStorageInternal()
308          */
309         bool requestStorage(int nframes);
310
311         /*! \brief
312          * Start storing data.
313          *
314          * \param[in] data    AbstractAnalysisData object containing this
315          *      storage.
316          * \param     modules Module manager for \p data.
317          * \exception std::bad_alloc if storage allocation fails.
318          *
319          * Typically called as \c startDataStorage(this, &moduleManager())
320          * from a member of \p data when the data is ready to be started.
321          * The storage object will take responsibility of calling all
322          * module notification methods in AnalysisDataModuleManager using
323          * \p modules.
324          *
325          * Lifetime of \p data and \p modules must exceed the lifetime of the
326          * storage object
327          * (typically, the storage object will be a member in \p data).
328          *
329          * Calls AnalysisDataModuleManager::notifyDataStart(), and throws any
330          * exceptions this method throws.
331          */
332         void startDataStorage(AbstractAnalysisData      *data,
333                               AnalysisDataModuleManager *modules);
334         /*! \brief
335          * Start storing data in parallel.
336          *
337          * \param[in] data    AbstractAnalysisData object containing this
338          *      storage.
339          * \param[in] options Parallelization options to use.
340          * \param     modules Module manager for \p data.
341          * \exception std::bad_alloc if storage allocation fails.
342          *
343          * Should be called instead of startDataStorage() if the data will be
344          * produced in parallel.  Works as startDataStorage(), but additionally
345          * initializes the storage and the attached modules to prepare for
346          * out-of-order data frames.
347          *
348          * Calls AnalysisDataModuleManager::notifyParallelDataStart(), and
349          * throws any exceptions this method throws.
350          */
351         void startParallelDataStorage(
352             AbstractAnalysisData              *data,
353             AnalysisDataModuleManager         *modules,
354             const AnalysisDataParallelOptions &options);
355         /*! \brief
356          * Starts storing a new frame.
357          *
358          * \param[in] header  Header for the new frame.
359          * \retval  Frame object corresponding to the started frame.
360          * \exception std::bad_alloc if storage reallocation fails
361          *      (only possible if storage of all frames has been requested).
362          * \exception APIError if frame is too far in the future.
363          *
364          * The returned object will be valid until the corresponding
365          * finishFrame() call.
366          *
367          * Must be called exactly once for each frame index.
368          *
369          * Currently, the implementation only works if the new frame is not too
370          * far in the future:
371          * If \c i is the index of the last frame such that all frames from
372          * 0, ..., \c i have been finished, then \p header().index() should be
373          * at most \c 2*parallelizationFactor-1 larger than \c i, where
374          * parallelizationFactor is the parallelization factor passed to
375          * setParallelOptions().
376          * Throws APIError if this constraint is violated.
377          *
378          * Calls AnalysisDataModuleManager::notifyFrameStart() (in certain
379          * cases) and AnalysisDataModuleManager::notifyParallelFrameStart(),
380          * and throws any exceptions these methods throw.
381          */
382         AnalysisDataStorageFrame &startFrame(const AnalysisDataFrameHeader &header);
383         /*! \brief
384          * Convenience method to start storing a new frame.
385          *
386          * Identical to \c startFrame(AnalysisDataFrameHeader(index, x, dx));
387          */
388         AnalysisDataStorageFrame &startFrame(int index, real x, real dx);
389         /*! \brief
390          * Obtain a frame object for an in-progress frame.
391          *
392          * \param[in] index  Frame index.
393          * \retval  Frame object corresponding to \p index.
394          *
395          * startFrame() should have been called for the frame with index
396          * \p index, and finishFrame() should not yet have been called.
397          * Returns the same object as returned by the original startFrame()
398          * call for the same index.
399          *
400          * Does not throw.
401          */
402         AnalysisDataStorageFrame &currentFrame(int index);
403         /*! \brief
404          * Convenience method for finishing a data frame.
405          *
406          * \param[in] index  Frame index.
407          *
408          * Identical to \c currentFrame(index).finishFrame().
409          *
410          * \see AnalysisDataStorageFrame::finishFrame()
411          */
412         void finishFrame(int index);
413         /*! \brief
414          * Finishes storing data.
415          *
416          * Calls AnalysisDataModuleManager::notifyDataFinish(), and throws any
417          * exceptions this method throws.
418          */
419         void finishDataStorage();
420
421     private:
422         typedef internal::AnalysisDataStorageImpl Impl;
423
424         PrivateImplPointer<Impl> impl_;
425 };
426
427 } // namespace gmx
428
429 #endif