Apply clang-format-11
[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), valueCount_(valueCount), dataSetIndex_(dataSetIndex), firstColumn_(firstColumn)
294     {
295         GMX_ASSERT(valueOffset >= 0, "Negative value offsets are invalid");
296         GMX_ASSERT(valueCount >= 0, "Negative value counts are invalid");
297         GMX_ASSERT(dataSetIndex >= 0, "Negative data set indices are invalid");
298         GMX_ASSERT(firstColumn >= 0, "Negative column indices are invalid");
299     }
300
301     //! Returns the offset of the first value in the referenced value array.
302     int valueOffset() const { return valueOffset_; }
303     //! Returns the number of values in this point set.
304     int valueCount() const { return valueCount_; }
305     //! Returns the data set index for this point set.
306     int dataSetIndex() const { return dataSetIndex_; }
307     //! Returns the index of the first column in this point set.
308     int firstColumn() const { return firstColumn_; }
309
310 private:
311     int valueOffset_;
312     int valueCount_;
313     int dataSetIndex_;
314     int firstColumn_;
315 };
316
317 //! Shorthand for reference to an array of point set data objects.
318 typedef ArrayRef<const AnalysisDataPointSetInfo> AnalysisDataPointSetInfosRef;
319
320 //! \endcond
321
322
323 /*! \brief
324  * Value type wrapper for non-mutable access to a set of data column values.
325  *
326  * Default copy constructor and assignment operator are used and work as
327  * intended.
328  * Typically new objects of this type are only constructed internally by the
329  * library and in classes that are derived from AbstractAnalysisData.
330  *
331  * Methods in this class do not throw, but may contain asserts for incorrect
332  * usage.
333  *
334  * The design of the interfaces is such that all objects of this type should be
335  * valid, i.e., header().isValid() should always return true.
336  *
337  * Note that it is not possible to change the contents of an initialized
338  * object, except by assigning a new object to replace it completely.
339  *
340  * \inpublicapi
341  * \ingroup module_analysisdata
342  */
343 class AnalysisDataPointSetRef
344 {
345 public:
346     /*! \brief
347      * Constructs a point set reference from given values.
348      *
349      * \param[in] header       Header for the frame.
350      * \param[in] pointSetInfo Information about the point set.
351      * \param[in] values       Values for each column.
352      *
353      * The first element of the point set should be found from \p values
354      * using the offset in \p pointSetInfo.
355      */
356     AnalysisDataPointSetRef(const AnalysisDataFrameHeader&  header,
357                             const AnalysisDataPointSetInfo& pointSetInfo,
358                             const AnalysisDataValuesRef&    values);
359     /*! \brief
360      * Constructs a point set reference from given values.
361      *
362      * \param[in] header      Header for the frame.
363      * \param[in] values      Values for each column.
364      *
365      * The first element in \p values should correspond to the first
366      * column.
367      */
368     AnalysisDataPointSetRef(const AnalysisDataFrameHeader&        header,
369                             const std::vector<AnalysisDataValue>& values);
370     /*! \brief
371      * Constructs a point set reference to a subset of columns.
372      *
373      * \param[in] points      Point set to use as source.
374      * \param[in] firstColumn First column index to include.
375      * \param[in] columnCount Number of columns to include.
376      *
377      * Creates a point set that contains \p columnCount columns starting
378      * from \p firstColumn from \p points, or a subset if all requested
379      * columns are not present in \p points.  If the requested column range
380      * and the range in \p points do not intersect, the result has
381      * columnCount() == 0.
382      *
383      * \p firstColumn is relative to the whole data set, i.e., not relative
384      * to points.firstColumn().
385      *
386      * Mainly intended for internal use.
387      */
388     AnalysisDataPointSetRef(const AnalysisDataPointSetRef& points, int firstColumn, int columnCount);
389
390     /*! \brief
391      * Returns the frame header for the frame of this point set.
392      */
393     const AnalysisDataFrameHeader& header() const { return header_; }
394     //! \copydoc AnalysisDataFrameHeader::index()
395     int frameIndex() const { return header_.index(); }
396     //! \copydoc AnalysisDataFrameHeader::x()
397     real x() const { return header_.x(); }
398     //! \copydoc AnalysisDataFrameHeader::dx()
399     real dx() const { return header_.dx(); }
400     //! Returns zero-based index of the dataset that this set is part of.
401     int dataSetIndex() const { return dataSetIndex_; }
402     //! Returns zero-based index of the first column included in this set.
403     int firstColumn() const { return firstColumn_; }
404     //! Returns the number of columns included in this set.
405     int columnCount() const { return ssize(values()); }
406     //! Returns zero-based index of the last column included in this set (inclusive).
407     int lastColumn() const { return firstColumn_ + columnCount() - 1; }
408     /*! \brief
409      * Returns reference container for all values.
410      *
411      * First value in the returned container corresponds to firstColumn().
412      */
413     const AnalysisDataValuesRef& values() const { return values_; }
414     /*! \brief
415      * Returns data value for a column in this set.
416      *
417      * \param[in] i  Zero-based column index relative to firstColumn().
418      *     Should be >= 0 and < columnCount().
419      */
420     real y(int i) const
421     {
422         GMX_ASSERT(i >= 0 && i < columnCount(), "Out of range data access");
423         return values()[i].value();
424     }
425     /*! \brief
426      * Returns error estimate for a column in this set if applicable.
427      *
428      * \param[in] i  Zero-based column index relative to firstColumn().
429      *     Should be >= 0 and < columnCount().
430      *
431      * Currently, this method returns zero if the source data does not
432      * specify errors.
433      */
434     real dy(int i) const
435     {
436         GMX_ASSERT(i >= 0 && i < columnCount(), "Out of range data access");
437         return values()[i].error();
438     }
439     /*! \brief
440      * Returns whether a column is present in this set.
441      *
442      * \param[in] i  Zero-based column index relative to firstColumn().
443      *     Should be >= 0 and < columnCount().
444      *
445      * If present(i) returns false, it is depends on the source data
446      * whether y(i) and/or dy(i) are defined.
447      */
448     bool present(int i) const
449     {
450         GMX_ASSERT(i >= 0 && i < columnCount(), "Out of range data access");
451         return values()[i].isPresent();
452     }
453     /*! \brief
454      * Returns true if all points in this point set are present.
455      *
456      * That is, if present() would return true for all points.
457      */
458     bool allPresent() const;
459
460 private:
461     AnalysisDataFrameHeader header_;
462     int                     dataSetIndex_;
463     int                     firstColumn_;
464     AnalysisDataValuesRef   values_;
465 };
466
467
468 /*! \brief
469  * Value type wrapper for non-mutable access to a data frame.
470  *
471  * Default copy constructor and assignment operator are used and work as
472  * intended.
473  * Typically new objects of this type are only constructed internally by the
474  * library and in classes that are derived from AbstractAnalysisData.
475  *
476  * Methods in this class do not throw, but may contain asserts for incorrect
477  * usage.
478  *
479  * Note that it is not possible to change the contents of an initialized
480  * object, except by assigning a new object to replace it completely.
481  *
482  * \inpublicapi
483  * \ingroup module_analysisdata
484  */
485 class AnalysisDataFrameRef
486 {
487 public:
488     /*! \brief
489      * Constructs an invalid frame reference.
490      *
491      * Return values of other methods than isValid() are unspecified for
492      * the constructed object.
493      */
494     AnalysisDataFrameRef();
495     /*! \brief
496      * Constructs a frame reference from given values.
497      *
498      * \param[in] header      Header for the frame.
499      * \param[in] values      Values for each column.
500      * \param[in] pointSets   Point set data.
501      */
502     AnalysisDataFrameRef(const AnalysisDataFrameHeader&      header,
503                          const AnalysisDataValuesRef&        values,
504                          const AnalysisDataPointSetInfosRef& pointSets);
505     /*! \brief
506      * Constructs a frame reference from given values.
507      *
508      * \param[in] header      Header for the frame.
509      * \param[in] values      Values for each column.
510      * \param[in] pointSets   Point set data.
511      */
512     AnalysisDataFrameRef(const AnalysisDataFrameHeader&               header,
513                          const std::vector<AnalysisDataValue>&        values,
514                          const std::vector<AnalysisDataPointSetInfo>& pointSets);
515     /*! \brief
516      * Constructs a frame reference to a subset of columns.
517      *
518      * \param[in] frame       Frame to use as source.
519      * \param[in] firstColumn First column index to include.
520      * \param[in] columnCount Number of columns to include.
521      *
522      * Creates a frame reference that contains \p columnCount columns
523      * starting from \p firstColumn from \p frame, or a subset if all
524      * requested columns are not present in \p frame.
525      *
526      * Mainly intended for internal use.
527      */
528     AnalysisDataFrameRef(const AnalysisDataFrameRef& frame, int firstColumn, int columnCount);
529
530     /*! \brief
531      * Returns whether the object refers to a valid frame.
532      *
533      * If returns false, return values of other methods are not specified.
534      */
535     bool isValid() const { return header().isValid(); }
536     //! Returns the header for this frame.
537     const AnalysisDataFrameHeader& header() const { return header_; }
538     //! \copydoc AnalysisDataFrameHeader::index()
539     int frameIndex() const { return header().index(); }
540     //! \copydoc AnalysisDataFrameHeader::x()
541     real x() const { return header().x(); }
542     //! \copydoc AnalysisDataFrameHeader::dx()
543     real dx() const { return header().dx(); }
544     /*! \brief
545      * Returns the number of point sets for this frame.
546      *
547      * Returns zero for an invalid frame.
548      */
549     int pointSetCount() const { return ssize(pointSets_); }
550     /*! \brief
551      * Returns point set reference for a given point set.
552      *
553      * Should not be called for invalid frames.
554      */
555     AnalysisDataPointSetRef pointSet(int index) const
556     {
557         GMX_ASSERT(isValid(), "Invalid data frame accessed");
558         GMX_ASSERT(index >= 0 && index < pointSetCount(), "Out of range data access");
559         return AnalysisDataPointSetRef(header_, pointSets_[index], values_);
560     }
561     /*! \brief
562      * Convenience method for accessing a column value in simple data.
563      *
564      * \copydetails AnalysisDataPointSetRef::y()
565      */
566     real y(int i) const { return singleColumnValue(i).value(); }
567     /*! \brief
568      * Convenience method for accessing error for a column value in simple
569      * data.
570      *
571      * \copydetails AnalysisDataPointSetRef::dy()
572      */
573     real dy(int i) const { return singleColumnValue(i).error(); }
574     /*! \brief
575      * Convenience method for accessing present status for a column in
576      * simple data.
577      *
578      * \copydetails AnalysisDataPointSetRef::present()
579      */
580     bool present(int i) const { return singleColumnValue(i).isPresent(); }
581     /*! \brief
582      * Returns true if all points in this frame are present.
583      */
584     bool allPresent() const;
585
586 private:
587     //! Helper method for accessing single columns in simple data.
588     const AnalysisDataValue& singleColumnValue(int i) const
589     {
590         GMX_ASSERT(isValid(), "Invalid data frame accessed");
591         GMX_ASSERT(pointSets_.size() == 1U && pointSets_[0].firstColumn() == 0,
592                    "Convenience method not available for multiple point sets");
593         GMX_ASSERT(i >= 0 && i < ssize(values_), "Out of range data access");
594         return values_[i];
595     }
596
597     AnalysisDataFrameHeader      header_;
598     AnalysisDataValuesRef        values_;
599     AnalysisDataPointSetInfosRef pointSets_;
600 };
601
602 } // namespace gmx
603
604 #endif