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