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