2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2017,2019, 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.
37 * Declares classes for accessing data frame information.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
41 * \ingroup module_analysisdata
43 #ifndef GMX_ANALYSISDATA_DATAFRAME_H
44 #define GMX_ANALYSISDATA_DATAFRAME_H
48 #include "gromacs/utility/arrayref.h"
49 #include "gromacs/utility/flags.h"
50 #include "gromacs/utility/gmxassert.h"
51 #include "gromacs/utility/real.h"
57 * Value type for representing a single value in analysis data objects.
59 * Default copy constructor and assignment operator are used and work as
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
69 * \ingroup module_analysisdata
71 class AnalysisDataValue
75 * Constructs an unset value.
77 AnalysisDataValue() : value_(0.0), error_(0.0) {}
79 * Constructs a value object with the given value.
81 * The constructed object is marked as set and present.
83 explicit AnalysisDataValue(real value
)
84 : value_(value
), error_(0.0)
87 flags_
.set(efPresent
);
91 * Direct access to the value.
93 * Assigning a value to this does not mark the value as set; setValue()
94 * must be used for this.
96 real
&value() { return value_
; }
98 * Direct access to the error estimate.
100 * Assigning a value to this does not mark the error estimate as set;
101 * setValue() must be used for this.
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_
; }
109 * Returns whether this value has been set.
111 * If this method returns false, the return value of value() and
112 * error() are undefined.
114 bool isSet() const { return flags_
.test(efSet
); }
116 * Returns whether the error estimate for this value has been set.
118 * If this method returns false, but isSet() returns true, error()
121 bool hasError() const { return flags_
.test(efErrorSet
); }
123 * Returns whether this value has been marked as present.
125 * If this method returns false, it is up to the source data to define
126 * whether isSet() may return true.
128 bool isPresent() const { return flags_
.test(efPresent
); }
130 //! Clears and unsets this value.
133 *this = AnalysisDataValue();
136 void setValue(real value
, bool bPresent
= true)
140 flags_
.set(efPresent
, bPresent
);
142 //! Sets this value and its error estimate.
143 void setValue(real value
, real error
, bool bPresent
= true)
148 flags_
.set(efErrorSet
);
149 flags_
.set(efPresent
, bPresent
);
151 //! Set only error estimate for this value.
152 void setError(real error
)
155 flags_
.set(efErrorSet
);
159 //! Possible flags for \a flags_.
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.
167 //! Value for this value.
169 //! Error estimate for this value, zero if not set.
171 //! Status flags for thise value.
172 FlagsTemplate
<Flag
> flags_
;
175 //! Shorthand for reference to an array of data values.
176 typedef ArrayRef
<const AnalysisDataValue
> AnalysisDataValuesRef
;
180 * Value type for storing frame-level information for analysis data.
182 * Default copy constructor and assignment operator are used and work as
184 * Typically new objects of this type are only constructed internally by the
185 * library and in classes that are derived from AbstractAnalysisData.
187 * Methods in this class do not throw, but may contain asserts for incorrect
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.
194 * \ingroup module_analysisdata
196 class AnalysisDataFrameHeader
200 * Constructs an invalid frame header.
202 * Return values of other methods than isValid() are unspecified for
203 * the constructed object.
205 AnalysisDataFrameHeader();
207 * Constructs a frame header from given values.
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.
213 AnalysisDataFrameHeader(int index
, real x
, real dx
);
216 * Returns whether the frame header corresponds to a valid frame.
218 * If returns false, return values of other methods are not specified.
225 * Returns zero-based index of the frame.
227 * The return value is >= 0 for valid frames.
228 * Should not be called for invalid frames.
232 GMX_ASSERT(isValid(), "Tried to access invalid frame header");
236 * Returns the x coordinate for the frame.
238 * Should not be called for invalid frames.
242 GMX_ASSERT(isValid(), "Tried to access invalid frame header");
246 * Returns error in the x coordinate for the frame (if applicable).
248 * All data do not provide error estimates.
249 * Typically returns zero in those cases.
251 * Should not be called for invalid frames.
255 GMX_ASSERT(isValid(), "Tried to access invalid frame header");
266 /*! \cond libinternal */
267 /*! \libinternal \brief
268 * Value type for internal indexing of point sets.
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
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
280 * Default copy constructor and assignment operator are used and work as
282 * Typically new objects of this type are only constructed internally by the
283 * library and in classes that are derived from AbstractAnalysisData.
285 * Methods in this class do not throw, but may contain asserts for incorrect
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.
292 * \ingroup module_analysisdata
294 class AnalysisDataPointSetInfo
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
)
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");
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_
; }
325 //! Shorthand for reference to an array of point set data objects.
326 typedef ArrayRef
<const AnalysisDataPointSetInfo
> AnalysisDataPointSetInfosRef
;
332 * Value type wrapper for non-mutable access to a set of data column values.
334 * Default copy constructor and assignment operator are used and work as
336 * Typically new objects of this type are only constructed internally by the
337 * library and in classes that are derived from AbstractAnalysisData.
339 * Methods in this class do not throw, but may contain asserts for incorrect
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.
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.
349 * \ingroup module_analysisdata
351 class AnalysisDataPointSetRef
355 * Constructs a point set reference from given values.
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.
361 * The first element of the point set should be found from \p values
362 * using the offset in \p pointSetInfo.
364 AnalysisDataPointSetRef(const AnalysisDataFrameHeader
&header
,
365 const AnalysisDataPointSetInfo
&pointSetInfo
,
366 const AnalysisDataValuesRef
&values
);
368 * Constructs a point set reference from given values.
370 * \param[in] header Header for the frame.
371 * \param[in] values Values for each column.
373 * The first element in \p values should correspond to the first
376 AnalysisDataPointSetRef(const AnalysisDataFrameHeader
&header
,
377 const std::vector
<AnalysisDataValue
> &values
);
379 * Constructs a point set reference to a subset of columns.
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.
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.
391 * \p firstColumn is relative to the whole data set, i.e., not relative
392 * to points.firstColumn().
394 * Mainly intended for internal use.
396 AnalysisDataPointSetRef(const AnalysisDataPointSetRef
&points
,
397 int firstColumn
, int columnCount
);
400 * Returns the frame header for the frame of this point set.
402 const AnalysisDataFrameHeader
&header() const
406 //! \copydoc AnalysisDataFrameHeader::index()
407 int frameIndex() const
409 return header_
.index();
411 //! \copydoc AnalysisDataFrameHeader::x()
416 //! \copydoc AnalysisDataFrameHeader::dx()
421 //! Returns zero-based index of the dataset that this set is part of.
422 int dataSetIndex() const
424 return dataSetIndex_
;
426 //! Returns zero-based index of the first column included in this set.
427 int firstColumn() const
431 //! Returns the number of columns included in this set.
432 int columnCount() const
434 return ssize(values());
436 //! Returns zero-based index of the last column included in this set (inclusive).
437 int lastColumn() const
439 return firstColumn_
+ columnCount() - 1;
442 * Returns reference container for all values.
444 * First value in the returned container corresponds to firstColumn().
446 const AnalysisDataValuesRef
&values() const
451 * Returns data value for a column in this set.
453 * \param[in] i Zero-based column index relative to firstColumn().
454 * Should be >= 0 and < columnCount().
458 GMX_ASSERT(i
>= 0 && i
< columnCount(), "Out of range data access");
459 return values()[i
].value();
462 * Returns error estimate for a column in this set if applicable.
464 * \param[in] i Zero-based column index relative to firstColumn().
465 * Should be >= 0 and < columnCount().
467 * Currently, this method returns zero if the source data does not
472 GMX_ASSERT(i
>= 0 && i
< columnCount(), "Out of range data access");
473 return values()[i
].error();
476 * Returns whether a column is present in this set.
478 * \param[in] i Zero-based column index relative to firstColumn().
479 * Should be >= 0 and < columnCount().
481 * If present(i) returns false, it is depends on the source data
482 * whether y(i) and/or dy(i) are defined.
484 bool present(int i
) const
486 GMX_ASSERT(i
>= 0 && i
< columnCount(), "Out of range data access");
487 return values()[i
].isPresent();
490 * Returns true if all points in this point set are present.
492 * That is, if present() would return true for all points.
494 bool allPresent() const;
497 AnalysisDataFrameHeader header_
;
500 AnalysisDataValuesRef values_
;
505 * Value type wrapper for non-mutable access to a data frame.
507 * Default copy constructor and assignment operator are used and work as
509 * Typically new objects of this type are only constructed internally by the
510 * library and in classes that are derived from AbstractAnalysisData.
512 * Methods in this class do not throw, but may contain asserts for incorrect
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.
519 * \ingroup module_analysisdata
521 class AnalysisDataFrameRef
525 * Constructs an invalid frame reference.
527 * Return values of other methods than isValid() are unspecified for
528 * the constructed object.
530 AnalysisDataFrameRef();
532 * Constructs a frame reference from given values.
534 * \param[in] header Header for the frame.
535 * \param[in] values Values for each column.
536 * \param[in] pointSets Point set data.
538 AnalysisDataFrameRef(const AnalysisDataFrameHeader
&header
,
539 const AnalysisDataValuesRef
&values
,
540 const AnalysisDataPointSetInfosRef
&pointSets
);
542 * Constructs a frame reference from given values.
544 * \param[in] header Header for the frame.
545 * \param[in] values Values for each column.
546 * \param[in] pointSets Point set data.
548 AnalysisDataFrameRef(const AnalysisDataFrameHeader
&header
,
549 const std::vector
<AnalysisDataValue
> &values
,
550 const std::vector
<AnalysisDataPointSetInfo
> &pointSets
);
552 * Constructs a frame reference to a subset of columns.
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.
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.
562 * Mainly intended for internal use.
564 AnalysisDataFrameRef(const AnalysisDataFrameRef
&frame
,
565 int firstColumn
, int columnCount
);
568 * Returns whether the object refers to a valid frame.
570 * If returns false, return values of other methods are not specified.
574 return header().isValid();
576 //! Returns the header for this frame.
577 const AnalysisDataFrameHeader
&header() const
581 //! \copydoc AnalysisDataFrameHeader::index()
582 int frameIndex() const
584 return header().index();
586 //! \copydoc AnalysisDataFrameHeader::x()
591 //! \copydoc AnalysisDataFrameHeader::dx()
594 return header().dx();
597 * Returns the number of point sets for this frame.
599 * Returns zero for an invalid frame.
601 int pointSetCount() const
603 return ssize(pointSets_
);
606 * Returns point set reference for a given point set.
608 * Should not be called for invalid frames.
610 AnalysisDataPointSetRef
pointSet(int index
) const
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_
);
618 * Convenience method for accessing a column value in simple data.
620 * \copydetails AnalysisDataPointSetRef::y()
624 return singleColumnValue(i
).value();
627 * Convenience method for accessing error for a column value in simple
630 * \copydetails AnalysisDataPointSetRef::dy()
634 return singleColumnValue(i
).error();
637 * Convenience method for accessing present status for a column in
640 * \copydetails AnalysisDataPointSetRef::present()
642 bool present(int i
) const
644 return singleColumnValue(i
).isPresent();
647 * Returns true if all points in this frame are present.
649 bool allPresent() const;
652 //! Helper method for accessing single columns in simple data.
653 const AnalysisDataValue
&singleColumnValue(int i
) const
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
< ssize(values_
),
659 "Out of range data access");
663 AnalysisDataFrameHeader header_
;
664 AnalysisDataValuesRef values_
;
665 AnalysisDataPointSetInfosRef pointSets_
;