Split AnalysisDataStorageFrame into two.
[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 AnalysisDataParallelOptions;
62
63 class AnalysisDataStorage;
64
65 namespace internal
66 {
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         //! Returns number of columns for the frame.
91         int columnCount() const { return values_.size(); }
92
93         /*! \brief
94          * Sets value for a column.
95          *
96          * \param[in] column  Zero-based column index.
97          * \param[in] value   Value to set for the column.
98          * \param[in] bPresent Present flag to set for the column.
99          *
100          * If called multiple times for a column (within one point set for
101          * multipoint data), old values are overwritten.
102          *
103          * Does not throw.
104          */
105         void setValue(int column, real value, bool bPresent = true)
106         {
107             GMX_ASSERT(column >= 0 && column < columnCount(),
108                        "Invalid column index");
109             values_[column].setValue(value, bPresent);
110         }
111         /*! \brief
112          * Sets value for a column.
113          *
114          * \param[in] column  Zero-based column index.
115          * \param[in] value   Value to set for the column.
116          * \param[in] error   Error estimate to set for the column.
117          * \param[in] bPresent Present flag to set for the column.
118          *
119          * If called multiple times for a column (within one point set for
120          * multipoint data), old values are overwritten.
121          *
122          * Does not throw.
123          */
124         void setValue(int column, real value, real error, bool bPresent = true)
125         {
126             GMX_ASSERT(column >= 0 && column < columnCount(),
127                        "Invalid column index");
128             values_[column].setValue(value, error, bPresent);
129         }
130         /*! \brief
131          * Access value for a column.
132          *
133          * \param[in] column  Zero-based column index.
134          *
135          * Should only be called after the column value has been set using
136          * setValue(); assigning a value to \c value(i) does not mark the
137          * column as set.
138          *
139          * Does not throw.
140          */
141         real &value(int column)
142         {
143             GMX_ASSERT(column >= 0 && column < columnCount(),
144                        "Invalid column index");
145             return values_[column].value();
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().
154          *
155          * Does not throw.
156          */
157         real value(int column) const
158         {
159             GMX_ASSERT(column >= 0 && column < columnCount(),
160                        "Invalid column index");
161             return values_[column].value();
162         }
163         /*! \brief
164          * Mark point set as finished for multipoint data.
165          *
166          * Must be called after each point set for multipoint data, including
167          * the last (i.e., no values must be set between the last call to this
168          * method and AnalysisDataStorage::finishFrame()).
169          * Must not be called for non-multipoint data.
170          *
171          * After this method has been called, all values appear as not set.
172          *
173          * Calls AbstractAnalysisData::notifyPointsAdd(), and throws any
174          * exception this method throws.
175          */
176         void finishPointSet();
177         /*! \brief
178          * Finish storing a frame.
179          *
180          * Must be called exactly once for each frame returned by startFrame(),
181          * after the corresponding call.
182          * The frame object must not be accessed after the call.
183          *
184          * Calls notification methods in AbstractAnalysisData, and throws any
185          * exceptions these methods throw.
186          */
187         void finishFrame();
188
189     private:
190
191         /*! \brief
192          * Create a new storage frame.
193          *
194          * \param[in] columnCount  Number of columns for the frame.
195          */
196         explicit AnalysisDataStorageFrame(int columnCount);
197
198         //! Clear all column values from the frame.
199         void clearValues();
200
201         //! Implementation data.
202         internal::AnalysisDataStorageFrameData *data_;
203         //! Values for the currently in-progress point set.
204         std::vector<AnalysisDataValue>          values_;
205
206         //! Needed for access to the constructor.
207         friend class AnalysisDataStorage;
208         //! Needed for managing the frame the object points to.
209         friend class internal::AnalysisDataStorageFrameData;
210
211         GMX_DISALLOW_COPY_AND_ASSIGN(AnalysisDataStorageFrame);
212 };
213
214 /*! \libinternal \brief
215  * Helper class that implements storage of data.
216  *
217  * This class implements a standard way of storing data to avoid implementing
218  * storage in each class derived from AbstractAnalysisData separately.
219  * To use this class in a class derived from AbstractAnalysisData, a member
220  * variable of this type should be declared and the data storage methods
221  * forwarded to tryGetDataFrame() and requestStorage() in that object.
222  * Storage properties should be set up, and then startDataStorage() called
223  * after calling AbstractAnalysisData::notifyDataStart().
224  * New frames can then be added using startFrame(), currentFrame() and
225  * finishFrame() methods.  These methods (and
226  * AnalysisDataStorageFrame::finishPointSet()) take the responsibility of
227  * calling AbstractAnalysisData::notifyFrameStart(),
228  * AbstractAnalysisData::notifyPointsAdd() and
229  * AbstractAnalysisData::notifyFrameFinish() appropriately.
230  *
231  * \todo
232  * Full support for multipoint data.
233  * Currently, multipoint data is only supported in serial pass-through mode
234  * without any storage.
235  *
236  * \todo
237  * Proper multi-threaded implementation.
238  *
239  * \inlibraryapi
240  * \ingroup module_analysisdata
241  */
242 class AnalysisDataStorage
243 {
244     public:
245         //! Constructs a storage object.
246         AnalysisDataStorage();
247         ~AnalysisDataStorage();
248
249         /*! \brief
250          * Sets whether the data will be multipoint.
251          *
252          * \exception APIError if storage has been requested.
253          *
254          * It is not mandatory to call this method (the multipoint property is
255          * automatically detected in startDataStorage()), but currently calling
256          * it will produce better diagnostic messages.
257          * When full support for multipoint data has been implemented, this
258          * method can be removed.
259          */
260         void setMultipoint(bool bMultipoint);
261         /*! \brief
262          * Set parallelization options for the storage.
263          *
264          * \param[in] opt  Parallization options to use.
265          *
266          * If this method is not called, the storage is set up for serial
267          * storage only.
268          *
269          * Does not throw.
270          */
271         void setParallelOptions(const AnalysisDataParallelOptions &opt);
272
273         /*! \brief
274          * Implements access to data frames.
275          *
276          * This method is designed such that calls to
277          * AbstractAnalysisData::tryGetDataFrameInternal() can be directly
278          * forwarded to this method.  See that method for more documentation.
279          *
280          * A valid reference for a frame will be returned after finishFrame()
281          * has been called for that frame.
282          *
283          * \see AbstractAnalysisData::tryGetDataFrameInternal()
284          */
285         AnalysisDataFrameRef tryGetDataFrame(int index) const;
286         /*! \brief
287          * Implements storage requests.
288          *
289          * This method is designed such that calls to
290          * AbstractAnalysisData::requestStorageInternal() can be directly
291          * forwarded to this method.  See that method for more documentation.
292          *
293          * Currently, multipoint data cannot be stored, but all other storage
294          * request will be honored.
295          *
296          * \see AbstractAnalysisData::requestStorageInternal()
297          */
298         bool requestStorage(int nframes);
299
300         /*! \brief
301          * Start storing data.
302          *
303          * \param  data  AbstractAnalysisData object containing this storage.
304          * \exception std::bad_alloc if storage allocation fails.
305          * \exception APIError if storage has been requested and \p data is
306          *      multipoint.
307          *
308          * Lifetime of \p data must exceed the lifetime of the storage object
309          * (typically, the storage object will be a member in \p data).
310          * The storage object will take responsibility of calling
311          * AbstractAnalysisData::notifyFrameStart(),
312          * AbstractAnalysisData::notifyPointsAdd() and
313          * AbstractAnalysisData::notifyFrameFinish() for \p data appropriately.
314          *
315          * AbstractAnalysisData::notifyDataStart() must have been called for
316          * \p data, because that may trigger storage requests from attached
317          * modules.
318          */
319         void startDataStorage(AbstractAnalysisData *data);
320         /*! \brief
321          * Starts storing a new frame.
322          *
323          * \param[in] header  Header for the new frame.
324          * \retval  Frame object corresponding to the started frame.
325          * \exception std::bad_alloc if storage reallocation fails
326          *      (only possible if storage of all frames has been requested).
327          * \exception APIError if frame is too far in the future.
328          *
329          * The returned object will be valid until the corresponding
330          * finishFrame() call.
331          *
332          * Must be called exactly once for each frame index.
333          *
334          * Currently, the implementation only works if the new frame is not too
335          * far in the future:
336          * If \c i is the index of the last frame such that all frames from
337          * 0, ..., \c i have been finished, then \p header().index() should be
338          * at most \c 2*parallelizationFactor-1 larger than \c i, where
339          * parallelizationFactor is the parallelization factor passed to
340          * setParallelOptions().
341          * Throws APIError if this constraint is violated.
342          *
343          * Calls AbstractAnalysisData::notifyDataStarted() in certain cases,
344          * and throws any exceptions this method throws.
345          */
346         AnalysisDataStorageFrame &startFrame(const AnalysisDataFrameHeader &header);
347         /*! \brief
348          * Convenience method to start storing a new frame.
349          *
350          * Identical to \c startFrame(AnalysisDataFrameHeader(index, x, dx));
351          */
352         AnalysisDataStorageFrame &startFrame(int index, real x, real dx);
353         /*! \brief
354          * Obtain a frame object for an in-progress frame.
355          *
356          * \param[in] index  Frame index.
357          * \retval  Frame object corresponding to \p index.
358          *
359          * startFrame() should have been called for the frame with index
360          * \p index, and finishFrame() should not yet have been called.
361          * Returns the same object as returned by the original startFrame()
362          * call for the same index.
363          *
364          * Does not throw.
365          */
366         AnalysisDataStorageFrame &currentFrame(int index);
367         /*! \brief
368          * Convenience method for finishing a data frame.
369          *
370          * \param[in] index  Frame index.
371          *
372          * Identical to \c currentFrame(index).finishFrame().
373          *
374          * \see AnalysisDataStorageFrame::finishFrame()
375          */
376         void finishFrame(int index);
377
378     private:
379         class Impl;
380
381         PrivateImplPointer<Impl> impl_;
382
383         //! Needed because the frame needs to access the implementation class.
384         friend class internal::AnalysisDataStorageFrameData;
385 };
386
387 } // namespace gmx
388
389 #endif