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