Merge branch 'release-4-6' into release-5-0
[alexxy/gromacs.git] / src / gromacs / analysisdata / dataframe.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  * 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 classes for accessing data frame information.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \inpublicapi
41  * \ingroup module_analysisdata
42  */
43 #ifndef GMX_ANALYSISDATA_DATAFRAME_H
44 #define GMX_ANALYSISDATA_DATAFRAME_H
45
46 #include <vector>
47
48 #include "../legacyheaders/types/simple.h"
49
50 #include "../utility/arrayref.h"
51 #include "../utility/flags.h"
52 #include "../utility/gmxassert.h"
53
54 namespace gmx
55 {
56
57 /*! \brief
58  * Value type for representing a single value in analysis data objects.
59  *
60  * Default copy constructor and assignment operator are used and work as
61  * intended.
62  *
63  * Methods in this class do not throw.
64  *
65  * Non-const methods are provided for use within the library only; currently
66  * it is not possible to access a non-const AnalysisDataValue through the
67  * public interface.
68  *
69  * \inpublicapi
70  * \ingroup module_analysisdata
71  */
72 class AnalysisDataValue
73 {
74     public:
75         /*! \brief
76          * Constructs an unset value.
77          */
78         AnalysisDataValue() : value_(0.0), error_(0.0) {}
79         /*! \brief
80          * Constructs a value object with the given value.
81          *
82          * The constructed object is marked as set and present.
83          */
84         explicit AnalysisDataValue(real value)
85             : value_(value), error_(0.0)
86         {
87             flags_.set(efSet);
88             flags_.set(efPresent);
89         }
90
91         /*! \brief
92          * Direct access to the value.
93          *
94          * Assigning a value to this does not mark the value as set; setValue()
95          * must be used for this.
96          */
97         real &value() { return value_; }
98         /*! \brief
99          * Direct access to the error estimate.
100          *
101          * Assigning a value to this does not mark the error estimate as set;
102          * setValue() must be used for this.
103          */
104         real &error() { return error_; }
105         //! Returns the value for this value.
106         real value() const { return value_; }
107         //! Returns the error estimate for this value, or zero if not set.
108         real error() const { return error_; }
109         /*! \brief
110          * Returns whether this value has been set.
111          *
112          * If this method returns false, the return value of value() and
113          * error() are undefined.
114          */
115         bool isSet() const { return flags_.test(efSet); }
116         /*! \brief
117          * Returns whether the error estimate for this value has been set.
118          *
119          * If this method returns false, but isSet() returns true, error()
120          * returns zero.
121          */
122         bool hasError() const { return flags_.test(efErrorSet); }
123         /*! \brief
124          * Returns whether this value has been marked as present.
125          *
126          * If this method returns false, it is up to the source data to define
127          * whether isSet() may return true.
128          */
129         bool isPresent() const { return flags_.test(efPresent); }
130
131         //! Clears and unsets this value.
132         void clear()
133         {
134             *this = AnalysisDataValue();
135         }
136         //! Sets this value.
137         void setValue(real value, bool bPresent = true)
138         {
139             value_ = value;
140             flags_.set(efSet);
141             flags_.set(efPresent, bPresent);
142         }
143         //! Sets this value and its error estimate.
144         void setValue(real value, real error, bool bPresent = true)
145         {
146             value_ = value;
147             error_ = error;
148             flags_.set(efSet);
149             flags_.set(efErrorSet);
150             flags_.set(efPresent, bPresent);
151         }
152         //! Set only error estimate for this value.
153         void setError(real error)
154         {
155             error_ = error;
156             flags_.set(efErrorSet);
157         }
158
159     private:
160         //! Possible flags for \a flags_.
161         enum Flag
162         {
163             efSet       = 1<<0, //!< Value has been set.
164             efErrorSet  = 1<<1, //!< Error estimate has been set.
165             efPresent   = 1<<2  //!< Value is set as present.
166         };
167
168         //! Value for this value.
169         real                    value_;
170         //! Error estimate for this value, zero if not set.
171         real                    error_;
172         //! Status flags for thise value.
173         FlagsTemplate<Flag>     flags_;
174 };
175
176 //! Shorthand for reference to an array of data values.
177 typedef ConstArrayRef<AnalysisDataValue> AnalysisDataValuesRef;
178
179
180 /*! \brief
181  * Value type for storing frame-level information for analysis data.
182  *
183  * Default copy constructor and assignment operator are used and work as
184  * intended.
185  * Typically new objects of this type are only constructed internally by the
186  * library and in classes that are derived from AbstractAnalysisData.
187  *
188  * Methods in this class do not throw, but may contain asserts for incorrect
189  * usage.
190  *
191  * Note that it is not possible to change the contents of an initialized
192  * object, except by assigning a new object to replace it completely.
193  *
194  * \inpublicapi
195  * \ingroup module_analysisdata
196  */
197 class AnalysisDataFrameHeader
198 {
199     public:
200         /*! \brief
201          * Constructs an invalid frame header.
202          *
203          * Return values of other methods than isValid() are unspecified for
204          * the constructed object.
205          */
206         AnalysisDataFrameHeader();
207         /*! \brief
208          * Constructs a frame header from given values.
209          *
210          * \param[in] index  Index of the frame. Must be >= 0.
211          * \param[in] x      x coordinate for the frame.
212          * \param[in] dx     Error estimate for x.
213          */
214         AnalysisDataFrameHeader(int index, real x, real dx);
215
216         /*! \brief
217          * Returns whether the frame header corresponds to a valid frame.
218          *
219          * If returns false, return values of other methods are not specified.
220          */
221         bool isValid() const
222         {
223             return index_ >= 0;
224         }
225         /*! \brief
226          * Returns zero-based index of the frame.
227          *
228          * The return value is >= 0 for valid frames.
229          * Should not be called for invalid frames.
230          */
231         int index() const
232         {
233             GMX_ASSERT(isValid(), "Tried to access invalid frame header");
234             return index_;
235         }
236         /*! \brief
237          * Returns the x coordinate for the frame.
238          *
239          * Should not be called for invalid frames.
240          */
241         real x() const
242         {
243             GMX_ASSERT(isValid(), "Tried to access invalid frame header");
244             return x_;
245         }
246         /*! \brief
247          * Returns error in the x coordinate for the frame (if applicable).
248          *
249          * All data do not provide error estimates.
250          * Typically returns zero in those cases.
251          *
252          * Should not be called for invalid frames.
253          */
254         real dx() const
255         {
256             GMX_ASSERT(isValid(), "Tried to access invalid frame header");
257             return dx_;
258         }
259
260     private:
261         int                     index_;
262         real                    x_;
263         real                    dx_;
264 };
265
266
267 /*! \cond libinternal */
268 /*! \libinternal \brief
269  * Value type for internal indexing of point sets.
270  *
271  * This class contains the necessary data to split an array of
272  * AnalysisDataValue objects into point sets.  It is always specified in the
273  * context of an array of AnalysisDataValues: the point set specified by this
274  * class contains valueCount() values, starting from the array index
275  * valueOffset().
276  * The value at location valueOffset() corresponds to column firstColumn().
277  * It is not necessary for code using the analysis data framework to know of
278  * this class, but it is declared in a public header to allow using it in other
279  * types.
280  *
281  * Default copy constructor and assignment operator are used and work as
282  * intended.
283  * Typically new objects of this type are only constructed internally by the
284  * library and in classes that are derived from AbstractAnalysisData.
285  *
286  * Methods in this class do not throw, but may contain asserts for incorrect
287  * usage.
288  *
289  * Note that it is not possible to change the contents of an initialized
290  * object, except by assigning a new object to replace it completely.
291  *
292  * \inlibraryapi
293  * \ingroup module_analysisdata
294  */
295 class AnalysisDataPointSetInfo
296 {
297     public:
298         //! Construct point set data object with the given values.
299         AnalysisDataPointSetInfo(int valueOffset, int valueCount,
300                                  int dataSetIndex, int firstColumn)
301             : valueOffset_(valueOffset), valueCount_(valueCount),
302               dataSetIndex_(dataSetIndex), firstColumn_(firstColumn)
303         {
304             GMX_ASSERT(valueOffset  >= 0, "Negative value offsets are invalid");
305             GMX_ASSERT(valueCount   >= 0, "Negative value counts are invalid");
306             GMX_ASSERT(dataSetIndex >= 0, "Negative data set indices are invalid");
307             GMX_ASSERT(firstColumn  >= 0, "Negative column indices are invalid");
308         }
309
310         //! Returns the offset of the first value in the referenced value array.
311         int valueOffset() const { return valueOffset_; }
312         //! Returns the number of values in this point set.
313         int valueCount() const { return valueCount_; }
314         //! Returns the data set index for this point set.
315         int dataSetIndex() const { return dataSetIndex_; }
316         //! Returns the index of the first column in this point set.
317         int firstColumn() const { return firstColumn_; }
318
319     private:
320         int                     valueOffset_;
321         int                     valueCount_;
322         int                     dataSetIndex_;
323         int                     firstColumn_;
324 };
325
326 //! Shorthand for reference to an array of point set data objects.
327 typedef ConstArrayRef<AnalysisDataPointSetInfo> AnalysisDataPointSetInfosRef;
328
329 //! \endcond
330
331
332 /*! \brief
333  * Value type wrapper for non-mutable access to a set of data column values.
334  *
335  * Default copy constructor and assignment operator are used and work as
336  * intended.
337  * Typically new objects of this type are only constructed internally by the
338  * library and in classes that are derived from AbstractAnalysisData.
339  *
340  * Methods in this class do not throw, but may contain asserts for incorrect
341  * usage.
342  *
343  * The design of the interfaces is such that all objects of this type should be
344  * valid, i.e., header().isValid() should always return true.
345  *
346  * Note that it is not possible to change the contents of an initialized
347  * object, except by assigning a new object to replace it completely.
348  *
349  * \inpublicapi
350  * \ingroup module_analysisdata
351  */
352 class AnalysisDataPointSetRef
353 {
354     public:
355         /*! \brief
356          * Constructs a point set reference from given values.
357          *
358          * \param[in] header       Header for the frame.
359          * \param[in] pointSetInfo Information about the point set.
360          * \param[in] values       Values for each column.
361          *
362          * The first element of the point set should be found from \p values
363          * using the offset in \p pointSetInfo.
364          */
365         AnalysisDataPointSetRef(const AnalysisDataFrameHeader  &header,
366                                 const AnalysisDataPointSetInfo &pointSetInfo,
367                                 const AnalysisDataValuesRef    &values);
368         /*! \brief
369          * Constructs a point set reference from given values.
370          *
371          * \param[in] header      Header for the frame.
372          * \param[in] values      Values for each column.
373          *
374          * The first element in \p values should correspond to the first
375          * column.
376          */
377         AnalysisDataPointSetRef(const AnalysisDataFrameHeader        &header,
378                                 const std::vector<AnalysisDataValue> &values);
379         /*! \brief
380          * Constructs a point set reference to a subset of columns.
381          *
382          * \param[in] points      Point set to use as source.
383          * \param[in] firstColumn First column index to include.
384          * \param[in] columnCount Number of columns to include.
385          *
386          * Creates a point set that contains \p columnCount columns starting
387          * from \p firstColumn from \p points, or a subset if all requested
388          * columns are not present in \p points.  If the requested column range
389          * and the range in \p points do not intersect, the result has
390          * columnCount() == 0.
391          *
392          * \p firstColumn is relative to the whole data set, i.e., not relative
393          * to points.firstColumn().
394          *
395          * Mainly intended for internal use.
396          */
397         AnalysisDataPointSetRef(const AnalysisDataPointSetRef &points,
398                                 int firstColumn, int columnCount);
399
400         /*! \brief
401          * Returns the frame header for the frame of this point set.
402          */
403         const AnalysisDataFrameHeader &header() const
404         {
405             return header_;
406         }
407         //! \copydoc AnalysisDataFrameHeader::index()
408         int frameIndex() const
409         {
410             return header_.index();
411         }
412         //! \copydoc AnalysisDataFrameHeader::x()
413         real x() const
414         {
415             return header_.x();
416         }
417         //! \copydoc AnalysisDataFrameHeader::dx()
418         real dx() const
419         {
420             return header_.dx();
421         }
422         //! Returns zero-based index of the dataset that this set is part of.
423         int dataSetIndex() const
424         {
425             return dataSetIndex_;
426         }
427         //! Returns zero-based index of the first column included in this set.
428         int firstColumn() const
429         {
430             return firstColumn_;
431         }
432         //! Returns the number of columns included in this set.
433         int columnCount() const
434         {
435             return values().size();
436         }
437         //! Returns zero-based index of the last column included in this set (inclusive).
438         int lastColumn() const
439         {
440             return firstColumn_ + columnCount() - 1;
441         }
442         /*! \brief
443          * Returns reference container for all values.
444          *
445          * First value in the returned container corresponds to firstColumn().
446          */
447         const AnalysisDataValuesRef &values() const
448         {
449             return values_;
450         }
451         /*! \brief
452          * Returns data value for a column in this set.
453          *
454          * \param[in] i  Zero-based column index relative to firstColumn().
455          *     Should be >= 0 and < columnCount().
456          */
457         real y(int i) const
458         {
459             GMX_ASSERT(i >= 0 && i < columnCount(), "Out of range data access");
460             return values()[i].value();
461         }
462         /*! \brief
463          * Returns error estimate for a column in this set if applicable.
464          *
465          * \param[in] i  Zero-based column index relative to firstColumn().
466          *     Should be >= 0 and < columnCount().
467          *
468          * Currently, this method returns zero if the source data does not
469          * specify errors.
470          */
471         real dy(int i) const
472         {
473             GMX_ASSERT(i >= 0 && i < columnCount(), "Out of range data access");
474             return values()[i].error();
475         }
476         /*! \brief
477          * Returns whether a column is present in this set.
478          *
479          * \param[in] i  Zero-based column index relative to firstColumn().
480          *     Should be >= 0 and < columnCount().
481          *
482          * If present(i) returns false, it is depends on the source data
483          * whether y(i) and/or dy(i) are defined.
484          */
485         bool present(int i) const
486         {
487             GMX_ASSERT(i >= 0 && i < columnCount(), "Out of range data access");
488             return values()[i].isPresent();
489         }
490         /*! \brief
491          * Returns true if all points in this point set are present.
492          *
493          * That is, if present() would return true for all points.
494          */
495         bool allPresent() const;
496
497     private:
498         AnalysisDataFrameHeader header_;
499         int                     dataSetIndex_;
500         int                     firstColumn_;
501         AnalysisDataValuesRef   values_;
502 };
503
504
505 /*! \brief
506  * Value type wrapper for non-mutable access to a data frame.
507  *
508  * Default copy constructor and assignment operator are used and work as
509  * intended.
510  * Typically new objects of this type are only constructed internally by the
511  * library and in classes that are derived from AbstractAnalysisData.
512  *
513  * Methods in this class do not throw, but may contain asserts for incorrect
514  * usage.
515  *
516  * Note that it is not possible to change the contents of an initialized
517  * object, except by assigning a new object to replace it completely.
518  *
519  * \inpublicapi
520  * \ingroup module_analysisdata
521  */
522 class AnalysisDataFrameRef
523 {
524     public:
525         /*! \brief
526          * Constructs an invalid frame reference.
527          *
528          * Return values of other methods than isValid() are unspecified for
529          * the constructed object.
530          */
531         AnalysisDataFrameRef();
532         /*! \brief
533          * Constructs a frame reference from given values.
534          *
535          * \param[in] header      Header for the frame.
536          * \param[in] values      Values for each column.
537          * \param[in] pointSets   Point set data.
538          */
539         AnalysisDataFrameRef(const AnalysisDataFrameHeader      &header,
540                              const AnalysisDataValuesRef        &values,
541                              const AnalysisDataPointSetInfosRef &pointSets);
542         /*! \brief
543          * Constructs a frame reference from given values.
544          *
545          * \param[in] header      Header for the frame.
546          * \param[in] values      Values for each column.
547          * \param[in] pointSets   Point set data.
548          */
549         AnalysisDataFrameRef(const AnalysisDataFrameHeader               &header,
550                              const std::vector<AnalysisDataValue>        &values,
551                              const std::vector<AnalysisDataPointSetInfo> &pointSets);
552         /*! \brief
553          * Constructs a frame reference to a subset of columns.
554          *
555          * \param[in] frame       Frame to use as source.
556          * \param[in] firstColumn First column index to include.
557          * \param[in] columnCount Number of columns to include.
558          *
559          * Creates a frame reference that contains \p columnCount columns
560          * starting from \p firstColumn from \p frame, or a subset if all
561          * requested columns are not present in \p frame.
562          *
563          * Mainly intended for internal use.
564          */
565         AnalysisDataFrameRef(const AnalysisDataFrameRef &frame,
566                              int firstColumn, int columnCount);
567
568         /*! \brief
569          * Returns whether the object refers to a valid frame.
570          *
571          * If returns false, return values of other methods are not specified.
572          */
573         bool isValid() const
574         {
575             return header().isValid();
576         }
577         //! Returns the header for this frame.
578         const AnalysisDataFrameHeader &header() const
579         {
580             return header_;
581         }
582         //! \copydoc AnalysisDataFrameHeader::index()
583         int frameIndex() const
584         {
585             return header().index();
586         }
587         //! \copydoc AnalysisDataFrameHeader::x()
588         real x() const
589         {
590             return header().x();
591         }
592         //! \copydoc AnalysisDataFrameHeader::dx()
593         real dx() const
594         {
595             return header().dx();
596         }
597         /*! \brief
598          * Returns the number of point sets for this frame.
599          *
600          * Returns zero for an invalid frame.
601          */
602         int pointSetCount() const
603         {
604             return pointSets_.size();
605         }
606         /*! \brief
607          * Returns point set reference for a given point set.
608          *
609          * Should not be called for invalid frames.
610          */
611         AnalysisDataPointSetRef pointSet(int index) const
612         {
613             GMX_ASSERT(isValid(), "Invalid data frame accessed");
614             GMX_ASSERT(index >= 0 && index < pointSetCount(),
615                        "Out of range data access");
616             return AnalysisDataPointSetRef(header_, pointSets_[index], values_);
617         }
618         /*! \brief
619          * Convenience method for accessing a column value in simple data.
620          *
621          * \copydetails AnalysisDataPointSetRef::y()
622          */
623         real y(int i) const
624         {
625             return singleColumnValue(i).value();
626         }
627         /*! \brief
628          * Convenience method for accessing error for a column value in simple
629          * data.
630          *
631          * \copydetails AnalysisDataPointSetRef::dy()
632          */
633         real dy(int i) const
634         {
635             return singleColumnValue(i).error();
636         }
637         /*! \brief
638          * Convenience method for accessing present status for a column in
639          * simple data.
640          *
641          * \copydetails AnalysisDataPointSetRef::present()
642          */
643         bool present(int i) const
644         {
645             return singleColumnValue(i).isPresent();
646         }
647         /*! \brief
648          * Returns true if all points in this frame are present.
649          */
650         bool allPresent() const;
651
652     private:
653         //! Helper method for accessing single columns in simple data.
654         const AnalysisDataValue &singleColumnValue(int i) const
655         {
656             GMX_ASSERT(isValid(), "Invalid data frame accessed");
657             GMX_ASSERT(pointSets_.size() == 1U && pointSets_[0].firstColumn() == 0,
658                        "Convenience method not available for multiple point sets");
659             GMX_ASSERT(i >= 0 && i < static_cast<int>(values_.size()),
660                        "Out of range data access");
661             return values_[i];
662         }
663
664         AnalysisDataFrameHeader      header_;
665         AnalysisDataValuesRef        values_;
666         AnalysisDataPointSetInfosRef pointSets_;
667 };
668
669 } // namespace gmx
670
671 #endif