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