c04ae1807273ea6819abadcd652f6bd03207bfbe
[alexxy/gromacs.git] / src / gromacs / analysisdata / analysisdata.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, 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 /*! \file
36  * \brief
37  * Declares gmx::AnalysisData and gmx::AnalysisDataHandle.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \inpublicapi
41  * \ingroup module_analysisdata
42  */
43 #ifndef GMX_ANALYSISDATA_ANALYSISDATA_H
44 #define GMX_ANALYSISDATA_ANALYSISDATA_H
45
46 #include "gromacs/analysisdata/abstractdata.h"
47 #include "gromacs/utility/real.h"
48
49 namespace gmx
50 {
51
52 class AnalysisDataHandle;
53 class AnalysisDataParallelOptions;
54
55 /*! \brief
56  * Parallelizable data container for raw data.
57  *
58  * This is the main class used to implement parallelizable data processing in
59  * analysis tools.  It is used by first creating an object and setting its
60  * properties using setDataSetCount(), setColumnCount() and setMultipoint(),
61  * and attaching necessary modules using addModule() etc.  Then one or more
62  * AnalysisDataHandle objects can be created using startData().  Each data
63  * handle can then be independently used to provide data frames (each frame
64  * must be provided by a single handle, but different frames can be freely
65  * mixed between the handles).  The finishFrameSerial() method must be called
66  * in serial for each frame, after one of the handles has been used to provide
67  * the data for that frame.  When all data has been provided, the handles
68  * are destroyed using finishData() (or AnalysisDataHandle::finishData()).
69  *
70  * When used through the trajectory analysis framework, calls to startData(),
71  * finishFrameSerial(), and finishData() are handled by the framework.
72  *
73  * \todo
74  * Parallel implementation is not complete.
75  *
76  * \if internal
77  * Special note for MPI implementation: assuming that the initialization of
78  * data objects is identical in all processes, associating the data objects
79  * in different MPI processes should be possible without changes in the
80  * interface.
81  * Alternative, more robust implementation could get a unique ID as parameter
82  * to the constructor or a separate function, but would require all tools to
83  * provide it.  With the current registration mechanism in
84  * TrajectoryAnalysisModule, this should be straightforward.
85  * \endif
86  *
87  * \inpublicapi
88  * \ingroup module_analysisdata
89  */
90 class AnalysisData : public AbstractAnalysisData
91 {
92     public:
93         /*! \brief
94          * Creates an empty analysis data object.
95          *
96          * \throws std::bad_alloc if out of memory.
97          */
98         AnalysisData();
99         virtual ~AnalysisData();
100
101         /*! \brief
102          * Sets the number of data sets.
103          *
104          * \param[in] dataSetCount  Number of data sets (must be > 0).
105          * \throws    std::bad_alloc if out of memory.
106          * \throws    APIError if modules have been added that are not
107          *      compatible with the new data set count.
108          *
109          * Must not be called after startData() has been called.
110          * If not called, a single data set is assumed.
111          * If called multiple times, the last call takes effect.
112          */
113         void setDataSetCount(int dataSetCount);
114         /*! \brief
115          * Sets the number of columns in a data set.
116          *
117          * \param[in] dataSet      Zero-based data set index.
118          * \param[in] columnCount  Number of columns in the data (must be > 0).
119          * \throws    APIError if modules have been added that are not
120          *      compatible with the new column count.
121          *
122          * Must be called before startData() for each data set.
123          * Must not be called after startData() has been called.
124          * If called multiple times for a data set, the last call takes effect.
125          */
126         void setColumnCount(int dataSet, int columnCount);
127         /*! \brief
128          * Sets whether the data contains multiple points per column per frame.
129          *
130          * \param[in] bMultipoint  Whether the data will allow multiple points
131          *      per column within a single frame.
132          * \throws    APIError if modules have been added that are not
133          *      compatible with the new setting.
134          *
135          * If this method is not called, the data is not multipoint.
136          *
137          * Must not be called after startData() has been called.
138          *
139          * \see isMultipoint()
140          */
141         void setMultipoint(bool bMultipoint);
142
143         virtual int frameCount() const;
144
145         /*! \brief
146          * Creates a handle for adding data.
147          *
148          * \param[in]  opt     Options for setting how this handle will be
149          *     used.
150          * \returns The created handle.
151          * \throws  std::bad_alloc if out of memory.
152          * \throws  APIError if any attached data module is not compatible.
153          * \throws  unspecified  Any exception thrown by attached data modules
154          *      in AnalysisDataModuleInterface::dataStarted().
155          *
156          * The caller should retain the returned handle (or a copy of it), and
157          * pass it to finishData() after successfully adding all data.
158          * The caller should discard the returned handle if an error occurs;
159          * memory allocated for the handle will be freed when the AnalysisData
160          * object is destroyed.
161          *
162          * The \p opt options should be the same for all calls to this method,
163          * and the number of calls should match the parallelization factor
164          * defined in \p opt.
165          */
166         AnalysisDataHandle startData(const AnalysisDataParallelOptions &opt);
167         /*! \brief
168          * Performs in-order sequential processing for the next frame.
169          *
170          * \param[in]  frameIndex Index of the frame that has been finished.
171          * \throws  unspecified  Any exception thrown by attached data modules
172          *      in AnalysisDataModuleInterface::frameFinishedSerial().
173          *
174          * This method should be called sequentially for each frame, after data
175          * for that frame has been produced.  It is not necessary to call this
176          * method if there is no parallelism, i.e., if only a single data
177          * handle is created and the parallelization options provided at that
178          * time do not indicate parallelism.
179          */
180         void finishFrameSerial(int frameIndex);
181         /*! \brief
182          * Destroys a handle after all data has been added.
183          *
184          * \param[in]  handle  Handle to destroy.
185          * \throws  unspecified  Any exception thrown by attached data modules
186          *      in AnalysisDataModuleInterface::dataFinished().
187          *
188          * \p handle must have been obtained from startData() of this object.
189          * The order of the calls with respect to the corresponding startData()
190          * calls is not important.
191          *
192          * The \p handle (and any copies) are invalid after the call.
193          */
194         void finishData(AnalysisDataHandle handle);
195
196     private:
197         virtual AnalysisDataFrameRef tryGetDataFrameInternal(int index) const;
198         virtual bool requestStorageInternal(int nframes);
199
200         class Impl;
201
202         PrivateImplPointer<Impl> impl_;
203
204         friend class AnalysisDataHandle;
205 };
206
207 namespace internal
208 {
209 class AnalysisDataHandleImpl;
210 }   // namespace internal
211
212 /*! \brief
213  * Handle for inserting data into AnalysisData.
214  *
215  * This class provides an interface for adding data frames into an AnalysisData
216  * object.  After a handle is obtained from AnalysisData::startData(), new
217  * frames can be added using startFrame().  Then values for that frame are set
218  * using provided methods (see below), and finishFrame() is called.  After all
219  * frames have been added, finishData() (or AnalysisData::finishData()) must be
220  * called.
221  *
222  * For simple (non-multipoint) data, within a frame values can be set using
223  * selectDataSet(), setPoint() and setPoints().  Setting the same column in the
224  * same data set multiple times overrides previously set values.
225  * When the frame is finished, attached modules are notified.
226  *
227  * Multipoint data works otherwise similarly, but requires finishPointSet() to
228  * be called for each set of points for which the modules need to be notified.
229  * Each point set starts empty (after startFrame() or finishPointSet()), and
230  * values can be set using setPoint()/setPoints().
231  * A single point set can contain values only for a single data set, which must
232  * be selected with selectDataSet() before setting any values.
233  * finishPointSet() must also be called for the last point set just before
234  * finishFrame().
235  *
236  * This class works like a pointer type: copying and assignment is lightweight,
237  * and all copies work interchangeably, accessing the same internal handle.
238  * However, normally you should only keep one copy of a handle, i.e., treat
239  * this type as movable.
240  * Several handles created from the same AnalysisData object can exist
241  * concurrently, but must currently operate on separate frames.
242  *
243  * \inpublicapi
244  * \ingroup module_analysisdata
245  */
246 class AnalysisDataHandle
247 {
248     public:
249         /*! \brief
250          * Constructs an invalid data handle.
251          *
252          * This constructor is provided for convenience in cases where it is
253          * easiest to declare an AnalysisDataHandle without immediately
254          * assigning a value to it.  Any attempt to call methods without first
255          * assigning a value from AnalysisData::startData() to the handle
256          * causes an assert.
257          *
258          * Does not throw.
259          */
260         AnalysisDataHandle();
261
262         //! Returns whether this data handle is valid.
263         bool isValid() const { return impl_ != NULL; }
264
265         /*! \brief
266          * Start data for a new frame.
267          *
268          * \param[in] index  Zero-based index for the frame to start.
269          * \param[in] x      x value for the frame.
270          * \param[in] dx     Error in x for the frame if applicable.
271          *
272          * \throws    unspecified  Any exception thrown by attached data
273          *      modules in AnalysisDataModuleInterface::frameStarted().
274          *
275          * Each \p index value 0, 1, ..., N (where N is the total number of
276          * frames) should be started exactly once by exactly one handle of an
277          * AnalysisData object.  The frames may be started out of order, but
278          * currently the implementation places some limitations on how far
279          * the index can be in the future (as counted from the first frame that
280          * is not finished).
281          */
282         void startFrame(int index, real x, real dx = 0.0);
283         /*! \brief
284          * Selects a data set for subsequent setPoint()/setPoints() calls.
285          *
286          * \param[in] index  Zero-based data set index.
287          *
288          * After startFrame(), the first data set is always selected.
289          * The set value is remembered until the end of the current frame, also
290          * across finishPointSet() calls.
291          *
292          * Does not throw.
293          */
294         void selectDataSet(int index);
295         /*! \brief
296          * Set a value for a single column for the current frame.
297          *
298          * \param[in] column  Zero-based column index.
299          * \param[in] value   Value to set for the column.
300          * \param[in] bPresent Present flag to set for the column.
301          *
302          * If called multiple times for a column (within one point set for
303          * multipoint data), old values are overwritten.
304          *
305          * Does not throw.
306          */
307         void setPoint(int column, real value, bool bPresent = true);
308         /*! \brief
309          * Set a value and its error estimate for a single column for the
310          * current frame.
311          *
312          * \param[in] column  Zero-based column index.
313          * \param[in] value   Value to set for the column.
314          * \param[in] error   Error estimate to set for the column.
315          * \param[in] bPresent Present flag to set for the column.
316          *
317          * If called multiple times for a column (within one point set for
318          * multipoint data), old values are overwritten.
319          *
320          * Does not throw.
321          */
322         void setPoint(int column, real value, real error, bool bPresent = true);
323         /*! \brief
324          * Set values for consecutive columns for the current frame.
325          *
326          * \param[in] firstColumn  Zero-based column index.
327          * \param[in] count        Number of columns to set.
328          * \param[in] values       Value array of \p column items.
329          *
330          * Equivalent to calling setPoint(firstColumn + i, values[i]) for
331          * i from 0 to count.
332          *
333          * Does not throw.
334          */
335         void setPoints(int firstColumn, int count, const real *values);
336         /*! \brief
337          * Finish data for the current point set.
338          *
339          * \throws    APIError if any attached data module is not compatible.
340          * \throws    unspecified  Any exception thrown by attached data
341          *      modules in AnalysisDataModuleInterface::pointsAdded().
342          *
343          * Must be called after each point set for multipoint data, including
344          * the last (i.e., no values must be set between the last call to this
345          * method and AnalysisDataStorage::finishFrame()).
346          * Must not be called for non-multipoint data.
347          */
348         void finishPointSet();
349         /*! \brief
350          * Finish data for the current frame.
351          *
352          * \throws    APIError if any attached data module is not compatible.
353          * \throws    unspecified  Any exception thrown by attached data
354          *      modules in frame notification methods.
355          */
356         void finishFrame();
357         //! Calls AnalysisData::finishData() for this handle.
358         void finishData();
359
360     private:
361         /*! \brief
362          * Creates a new data handle associated with \p data.
363          *
364          * \param  impl Data to associate the handle with.
365          *
366          * The constructor is private because data handles should only be
367          * constructed through AnalysisData::startData().
368          *
369          * Does not throw.
370          */
371         explicit AnalysisDataHandle(internal::AnalysisDataHandleImpl *impl);
372
373         /*! \brief
374          * Pointer to the internal implementation class.
375          *
376          * The memory for this object is managed by the AnalysisData object,
377          * and AnalysisDataHandle simply provides a public interface for
378          * accessing the implementation.
379          */
380         internal::AnalysisDataHandleImpl *impl_;
381
382         /*! \brief
383          * Needed to access the non-public implementation.
384          */
385         friend class AnalysisData;
386 };
387
388 } // namespace gmx
389
390 #endif