Update instructions in containers.rst
[gromacs.git] / src / gromacs / analysisdata / dataframe.h
blobc034fd8ed27a030d85e12cc7e7afe5c74f13f566
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2017,2019,2020, 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.
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.
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.
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.
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.
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.
35 /*! \file
36 * \brief
37 * Declares classes for accessing data frame information.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \inpublicapi
41 * \ingroup module_analysisdata
43 #ifndef GMX_ANALYSISDATA_DATAFRAME_H
44 #define GMX_ANALYSISDATA_DATAFRAME_H
46 #include <vector>
48 #include "gromacs/utility/arrayref.h"
49 #include "gromacs/utility/flags.h"
50 #include "gromacs/utility/gmxassert.h"
51 #include "gromacs/utility/real.h"
53 namespace gmx
56 /*! \brief
57 * Value type for representing a single value in analysis data objects.
59 * Default copy constructor and assignment operator are used and work as
60 * intended.
62 * Methods in this class do not throw.
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.
68 * \inpublicapi
69 * \ingroup module_analysisdata
71 class AnalysisDataValue
73 public:
74 /*! \brief
75 * Constructs an unset value.
77 AnalysisDataValue() : value_(0.0), error_(0.0) {}
78 /*! \brief
79 * Constructs a value object with the given value.
81 * The constructed object is marked as set and present.
83 explicit AnalysisDataValue(real value) : value_(value), error_(0.0)
85 flags_.set(efSet);
86 flags_.set(efPresent);
89 /*! \brief
90 * Direct access to the value.
92 * Assigning a value to this does not mark the value as set; setValue()
93 * must be used for this.
95 real& value() { return value_; }
96 /*! \brief
97 * Direct access to the error estimate.
99 * Assigning a value to this does not mark the error estimate as set;
100 * setValue() must be used for this.
102 real& error() { return error_; }
103 //! Returns the value for this value.
104 real value() const { return value_; }
105 //! Returns the error estimate for this value, or zero if not set.
106 real error() const { return error_; }
107 /*! \brief
108 * Returns whether this value has been set.
110 * If this method returns false, the return value of value() and
111 * error() are undefined.
113 bool isSet() const { return flags_.test(efSet); }
114 /*! \brief
115 * Returns whether the error estimate for this value has been set.
117 * If this method returns false, but isSet() returns true, error()
118 * returns zero.
120 bool hasError() const { return flags_.test(efErrorSet); }
121 /*! \brief
122 * Returns whether this value has been marked as present.
124 * If this method returns false, it is up to the source data to define
125 * whether isSet() may return true.
127 bool isPresent() const { return flags_.test(efPresent); }
129 //! Clears and unsets this value.
130 void clear() { *this = AnalysisDataValue(); }
131 //! Sets this value.
132 void setValue(real value, bool bPresent = true)
134 value_ = value;
135 flags_.set(efSet);
136 flags_.set(efPresent, bPresent);
138 //! Sets this value and its error estimate.
139 void setValue(real value, real error, bool bPresent = true)
141 value_ = value;
142 error_ = error;
143 flags_.set(efSet);
144 flags_.set(efErrorSet);
145 flags_.set(efPresent, bPresent);
147 //! Set only error estimate for this value.
148 void setError(real error)
150 error_ = error;
151 flags_.set(efErrorSet);
154 private:
155 //! Possible flags for \a flags_.
156 enum Flag : uint64_t
158 efSet = 1 << 0, //!< Value has been set.
159 efErrorSet = 1 << 1, //!< Error estimate has been set.
160 efPresent = 1 << 2 //!< Value is set as present.
163 //! Value for this value.
164 real value_;
165 //! Error estimate for this value, zero if not set.
166 real error_;
167 //! Status flags for thise value.
168 FlagsTemplate<Flag> flags_;
171 //! Shorthand for reference to an array of data values.
172 typedef ArrayRef<const AnalysisDataValue> AnalysisDataValuesRef;
175 /*! \brief
176 * Value type for storing frame-level information for analysis data.
178 * Default copy constructor and assignment operator are used and work as
179 * intended.
180 * Typically new objects of this type are only constructed internally by the
181 * library and in classes that are derived from AbstractAnalysisData.
183 * Methods in this class do not throw, but may contain asserts for incorrect
184 * usage.
186 * Note that it is not possible to change the contents of an initialized
187 * object, except by assigning a new object to replace it completely.
189 * \inpublicapi
190 * \ingroup module_analysisdata
192 class AnalysisDataFrameHeader
194 public:
195 /*! \brief
196 * Constructs an invalid frame header.
198 * Return values of other methods than isValid() are unspecified for
199 * the constructed object.
201 AnalysisDataFrameHeader();
202 /*! \brief
203 * Constructs a frame header from given values.
205 * \param[in] index Index of the frame. Must be >= 0.
206 * \param[in] x x coordinate for the frame.
207 * \param[in] dx Error estimate for x.
209 AnalysisDataFrameHeader(int index, real x, real dx);
211 /*! \brief
212 * Returns whether the frame header corresponds to a valid frame.
214 * If returns false, return values of other methods are not specified.
216 bool isValid() const { return index_ >= 0; }
217 /*! \brief
218 * Returns zero-based index of the frame.
220 * The return value is >= 0 for valid frames.
221 * Should not be called for invalid frames.
223 int index() const
225 GMX_ASSERT(isValid(), "Tried to access invalid frame header");
226 return index_;
228 /*! \brief
229 * Returns the x coordinate for the frame.
231 * Should not be called for invalid frames.
233 real x() const
235 GMX_ASSERT(isValid(), "Tried to access invalid frame header");
236 return x_;
238 /*! \brief
239 * Returns error in the x coordinate for the frame (if applicable).
241 * All data do not provide error estimates.
242 * Typically returns zero in those cases.
244 * Should not be called for invalid frames.
246 real dx() const
248 GMX_ASSERT(isValid(), "Tried to access invalid frame header");
249 return dx_;
252 private:
253 int index_;
254 real x_;
255 real dx_;
259 /*! \cond libinternal */
260 /*! \libinternal \brief
261 * Value type for internal indexing of point sets.
263 * This class contains the necessary data to split an array of
264 * AnalysisDataValue objects into point sets. It is always specified in the
265 * context of an array of AnalysisDataValues: the point set specified by this
266 * class contains valueCount() values, starting from the array index
267 * valueOffset().
268 * The value at location valueOffset() corresponds to column firstColumn().
269 * It is not necessary for code using the analysis data framework to know of
270 * this class, but it is declared in a public header to allow using it in other
271 * types.
273 * Default copy constructor and assignment operator are used and work as
274 * intended.
275 * Typically new objects of this type are only constructed internally by the
276 * library and in classes that are derived from AbstractAnalysisData.
278 * Methods in this class do not throw, but may contain asserts for incorrect
279 * usage.
281 * Note that it is not possible to change the contents of an initialized
282 * object, except by assigning a new object to replace it completely.
284 * \inlibraryapi
285 * \ingroup module_analysisdata
287 class AnalysisDataPointSetInfo
289 public:
290 //! Construct point set data object with the given values.
291 AnalysisDataPointSetInfo(int valueOffset, int valueCount, int dataSetIndex, int firstColumn) :
292 valueOffset_(valueOffset),
293 valueCount_(valueCount),
294 dataSetIndex_(dataSetIndex),
295 firstColumn_(firstColumn)
297 GMX_ASSERT(valueOffset >= 0, "Negative value offsets are invalid");
298 GMX_ASSERT(valueCount >= 0, "Negative value counts are invalid");
299 GMX_ASSERT(dataSetIndex >= 0, "Negative data set indices are invalid");
300 GMX_ASSERT(firstColumn >= 0, "Negative column indices are invalid");
303 //! Returns the offset of the first value in the referenced value array.
304 int valueOffset() const { return valueOffset_; }
305 //! Returns the number of values in this point set.
306 int valueCount() const { return valueCount_; }
307 //! Returns the data set index for this point set.
308 int dataSetIndex() const { return dataSetIndex_; }
309 //! Returns the index of the first column in this point set.
310 int firstColumn() const { return firstColumn_; }
312 private:
313 int valueOffset_;
314 int valueCount_;
315 int dataSetIndex_;
316 int firstColumn_;
319 //! Shorthand for reference to an array of point set data objects.
320 typedef ArrayRef<const AnalysisDataPointSetInfo> AnalysisDataPointSetInfosRef;
322 //! \endcond
325 /*! \brief
326 * Value type wrapper for non-mutable access to a set of data column values.
328 * Default copy constructor and assignment operator are used and work as
329 * intended.
330 * Typically new objects of this type are only constructed internally by the
331 * library and in classes that are derived from AbstractAnalysisData.
333 * Methods in this class do not throw, but may contain asserts for incorrect
334 * usage.
336 * The design of the interfaces is such that all objects of this type should be
337 * valid, i.e., header().isValid() should always return true.
339 * Note that it is not possible to change the contents of an initialized
340 * object, except by assigning a new object to replace it completely.
342 * \inpublicapi
343 * \ingroup module_analysisdata
345 class AnalysisDataPointSetRef
347 public:
348 /*! \brief
349 * Constructs a point set reference from given values.
351 * \param[in] header Header for the frame.
352 * \param[in] pointSetInfo Information about the point set.
353 * \param[in] values Values for each column.
355 * The first element of the point set should be found from \p values
356 * using the offset in \p pointSetInfo.
358 AnalysisDataPointSetRef(const AnalysisDataFrameHeader& header,
359 const AnalysisDataPointSetInfo& pointSetInfo,
360 const AnalysisDataValuesRef& values);
361 /*! \brief
362 * Constructs a point set reference from given values.
364 * \param[in] header Header for the frame.
365 * \param[in] values Values for each column.
367 * The first element in \p values should correspond to the first
368 * column.
370 AnalysisDataPointSetRef(const AnalysisDataFrameHeader& header,
371 const std::vector<AnalysisDataValue>& values);
372 /*! \brief
373 * Constructs a point set reference to a subset of columns.
375 * \param[in] points Point set to use as source.
376 * \param[in] firstColumn First column index to include.
377 * \param[in] columnCount Number of columns to include.
379 * Creates a point set that contains \p columnCount columns starting
380 * from \p firstColumn from \p points, or a subset if all requested
381 * columns are not present in \p points. If the requested column range
382 * and the range in \p points do not intersect, the result has
383 * columnCount() == 0.
385 * \p firstColumn is relative to the whole data set, i.e., not relative
386 * to points.firstColumn().
388 * Mainly intended for internal use.
390 AnalysisDataPointSetRef(const AnalysisDataPointSetRef& points, int firstColumn, int columnCount);
392 /*! \brief
393 * Returns the frame header for the frame of this point set.
395 const AnalysisDataFrameHeader& header() const { return header_; }
396 //! \copydoc AnalysisDataFrameHeader::index()
397 int frameIndex() const { return header_.index(); }
398 //! \copydoc AnalysisDataFrameHeader::x()
399 real x() const { return header_.x(); }
400 //! \copydoc AnalysisDataFrameHeader::dx()
401 real dx() const { return header_.dx(); }
402 //! Returns zero-based index of the dataset that this set is part of.
403 int dataSetIndex() const { return dataSetIndex_; }
404 //! Returns zero-based index of the first column included in this set.
405 int firstColumn() const { return firstColumn_; }
406 //! Returns the number of columns included in this set.
407 int columnCount() const { return ssize(values()); }
408 //! Returns zero-based index of the last column included in this set (inclusive).
409 int lastColumn() const { return firstColumn_ + columnCount() - 1; }
410 /*! \brief
411 * Returns reference container for all values.
413 * First value in the returned container corresponds to firstColumn().
415 const AnalysisDataValuesRef& values() const { return values_; }
416 /*! \brief
417 * Returns data value for a column in this set.
419 * \param[in] i Zero-based column index relative to firstColumn().
420 * Should be >= 0 and < columnCount().
422 real y(int i) const
424 GMX_ASSERT(i >= 0 && i < columnCount(), "Out of range data access");
425 return values()[i].value();
427 /*! \brief
428 * Returns error estimate for a column in this set if applicable.
430 * \param[in] i Zero-based column index relative to firstColumn().
431 * Should be >= 0 and < columnCount().
433 * Currently, this method returns zero if the source data does not
434 * specify errors.
436 real dy(int i) const
438 GMX_ASSERT(i >= 0 && i < columnCount(), "Out of range data access");
439 return values()[i].error();
441 /*! \brief
442 * Returns whether a column is present in this set.
444 * \param[in] i Zero-based column index relative to firstColumn().
445 * Should be >= 0 and < columnCount().
447 * If present(i) returns false, it is depends on the source data
448 * whether y(i) and/or dy(i) are defined.
450 bool present(int i) const
452 GMX_ASSERT(i >= 0 && i < columnCount(), "Out of range data access");
453 return values()[i].isPresent();
455 /*! \brief
456 * Returns true if all points in this point set are present.
458 * That is, if present() would return true for all points.
460 bool allPresent() const;
462 private:
463 AnalysisDataFrameHeader header_;
464 int dataSetIndex_;
465 int firstColumn_;
466 AnalysisDataValuesRef values_;
470 /*! \brief
471 * Value type wrapper for non-mutable access to a data frame.
473 * Default copy constructor and assignment operator are used and work as
474 * intended.
475 * Typically new objects of this type are only constructed internally by the
476 * library and in classes that are derived from AbstractAnalysisData.
478 * Methods in this class do not throw, but may contain asserts for incorrect
479 * usage.
481 * Note that it is not possible to change the contents of an initialized
482 * object, except by assigning a new object to replace it completely.
484 * \inpublicapi
485 * \ingroup module_analysisdata
487 class AnalysisDataFrameRef
489 public:
490 /*! \brief
491 * Constructs an invalid frame reference.
493 * Return values of other methods than isValid() are unspecified for
494 * the constructed object.
496 AnalysisDataFrameRef();
497 /*! \brief
498 * Constructs a frame reference from given values.
500 * \param[in] header Header for the frame.
501 * \param[in] values Values for each column.
502 * \param[in] pointSets Point set data.
504 AnalysisDataFrameRef(const AnalysisDataFrameHeader& header,
505 const AnalysisDataValuesRef& values,
506 const AnalysisDataPointSetInfosRef& pointSets);
507 /*! \brief
508 * Constructs a frame reference from given values.
510 * \param[in] header Header for the frame.
511 * \param[in] values Values for each column.
512 * \param[in] pointSets Point set data.
514 AnalysisDataFrameRef(const AnalysisDataFrameHeader& header,
515 const std::vector<AnalysisDataValue>& values,
516 const std::vector<AnalysisDataPointSetInfo>& pointSets);
517 /*! \brief
518 * Constructs a frame reference to a subset of columns.
520 * \param[in] frame Frame to use as source.
521 * \param[in] firstColumn First column index to include.
522 * \param[in] columnCount Number of columns to include.
524 * Creates a frame reference that contains \p columnCount columns
525 * starting from \p firstColumn from \p frame, or a subset if all
526 * requested columns are not present in \p frame.
528 * Mainly intended for internal use.
530 AnalysisDataFrameRef(const AnalysisDataFrameRef& frame, int firstColumn, int columnCount);
532 /*! \brief
533 * Returns whether the object refers to a valid frame.
535 * If returns false, return values of other methods are not specified.
537 bool isValid() const { return header().isValid(); }
538 //! Returns the header for this frame.
539 const AnalysisDataFrameHeader& header() const { return header_; }
540 //! \copydoc AnalysisDataFrameHeader::index()
541 int frameIndex() const { return header().index(); }
542 //! \copydoc AnalysisDataFrameHeader::x()
543 real x() const { return header().x(); }
544 //! \copydoc AnalysisDataFrameHeader::dx()
545 real dx() const { return header().dx(); }
546 /*! \brief
547 * Returns the number of point sets for this frame.
549 * Returns zero for an invalid frame.
551 int pointSetCount() const { return ssize(pointSets_); }
552 /*! \brief
553 * Returns point set reference for a given point set.
555 * Should not be called for invalid frames.
557 AnalysisDataPointSetRef pointSet(int index) const
559 GMX_ASSERT(isValid(), "Invalid data frame accessed");
560 GMX_ASSERT(index >= 0 && index < pointSetCount(), "Out of range data access");
561 return AnalysisDataPointSetRef(header_, pointSets_[index], values_);
563 /*! \brief
564 * Convenience method for accessing a column value in simple data.
566 * \copydetails AnalysisDataPointSetRef::y()
568 real y(int i) const { return singleColumnValue(i).value(); }
569 /*! \brief
570 * Convenience method for accessing error for a column value in simple
571 * data.
573 * \copydetails AnalysisDataPointSetRef::dy()
575 real dy(int i) const { return singleColumnValue(i).error(); }
576 /*! \brief
577 * Convenience method for accessing present status for a column in
578 * simple data.
580 * \copydetails AnalysisDataPointSetRef::present()
582 bool present(int i) const { return singleColumnValue(i).isPresent(); }
583 /*! \brief
584 * Returns true if all points in this frame are present.
586 bool allPresent() const;
588 private:
589 //! Helper method for accessing single columns in simple data.
590 const AnalysisDataValue& singleColumnValue(int i) const
592 GMX_ASSERT(isValid(), "Invalid data frame accessed");
593 GMX_ASSERT(pointSets_.size() == 1U && pointSets_[0].firstColumn() == 0,
594 "Convenience method not available for multiple point sets");
595 GMX_ASSERT(i >= 0 && i < ssize(values_), "Out of range data access");
596 return values_[i];
599 AnalysisDataFrameHeader header_;
600 AnalysisDataValuesRef values_;
601 AnalysisDataPointSetInfosRef pointSets_;
604 } // namespace gmx
606 #endif