Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / analysisdata / dataframe.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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 /*! \internal \file
36  * \brief
37  * Implements classes in dataframe.h.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_analysisdata
41  */
42 #include "gmxpre.h"
43
44 #include "dataframe.h"
45
46 #include "gromacs/utility/gmxassert.h"
47
48 namespace gmx
49 {
50
51 /********************************************************************
52  * AnalysisDataFrameHeader
53  */
54
55 AnalysisDataFrameHeader::AnalysisDataFrameHeader() : index_(-1), x_(0.0), dx_(0.0) {}
56
57
58 AnalysisDataFrameHeader::AnalysisDataFrameHeader(int index, real x, real dx) :
59     index_(index),
60     x_(x),
61     dx_(dx)
62 {
63     GMX_ASSERT(index >= 0, "Invalid frame index");
64 }
65
66
67 /********************************************************************
68  * AnalysisDataPointSetRef
69  */
70
71 AnalysisDataPointSetRef::AnalysisDataPointSetRef(const AnalysisDataFrameHeader&  header,
72                                                  const AnalysisDataPointSetInfo& pointSetInfo,
73                                                  const AnalysisDataValuesRef&    values) :
74     header_(header),
75     dataSetIndex_(pointSetInfo.dataSetIndex()),
76     firstColumn_(pointSetInfo.firstColumn()),
77     values_(constArrayRefFromArray(&*values.begin() + pointSetInfo.valueOffset(), pointSetInfo.valueCount()))
78 {
79     GMX_ASSERT(header_.isValid(), "Invalid point set reference should not be constructed");
80 }
81
82
83 AnalysisDataPointSetRef::AnalysisDataPointSetRef(const AnalysisDataFrameHeader&        header,
84                                                  const std::vector<AnalysisDataValue>& values) :
85     header_(header),
86     dataSetIndex_(0),
87     firstColumn_(0),
88     values_(values)
89 {
90     GMX_ASSERT(header_.isValid(), "Invalid point set reference should not be constructed");
91 }
92
93
94 AnalysisDataPointSetRef::AnalysisDataPointSetRef(const AnalysisDataPointSetRef& points,
95                                                  int                            firstColumn,
96                                                  int                            columnCount) :
97     header_(points.header()),
98     dataSetIndex_(points.dataSetIndex()),
99     firstColumn_(0)
100 {
101     GMX_ASSERT(firstColumn >= 0, "Invalid first column");
102     GMX_ASSERT(columnCount >= 0, "Invalid column count");
103     if (points.lastColumn() < firstColumn || points.firstColumn() >= firstColumn + columnCount
104         || columnCount == 0)
105     {
106         return;
107     }
108     AnalysisDataValuesRef::const_iterator begin        = points.values().begin();
109     int                                   pointsOffset = firstColumn - points.firstColumn();
110     if (pointsOffset > 0)
111     {
112         // Offset pointer if the first column is not the first in points.
113         begin += pointsOffset;
114     }
115     else
116     {
117         // Take into account if first column is before the first in points.
118         firstColumn_ = -pointsOffset;
119         columnCount -= -pointsOffset;
120     }
121     // Decrease column count if there are not enough columns in points.
122     AnalysisDataValuesRef::const_iterator end = begin + columnCount;
123     if (pointsOffset + columnCount > points.columnCount())
124     {
125         end = points.values().end();
126     }
127     values_ = AnalysisDataValuesRef(begin, end);
128 }
129
130
131 bool AnalysisDataPointSetRef::allPresent() const
132 {
133     AnalysisDataValuesRef::const_iterator i;
134     for (i = values_.begin(); i != values_.end(); ++i)
135     {
136         if (!i->isPresent())
137         {
138             return false;
139         }
140     }
141     return true;
142 }
143
144
145 /********************************************************************
146  * AnalysisDataFrameRef
147  */
148
149 AnalysisDataFrameRef::AnalysisDataFrameRef() {}
150
151
152 AnalysisDataFrameRef::AnalysisDataFrameRef(const AnalysisDataFrameHeader&      header,
153                                            const AnalysisDataValuesRef&        values,
154                                            const AnalysisDataPointSetInfosRef& pointSets) :
155     header_(header),
156     values_(values),
157     pointSets_(pointSets)
158 {
159     GMX_ASSERT(!pointSets_.empty(), "There must always be a point set");
160 }
161
162
163 AnalysisDataFrameRef::AnalysisDataFrameRef(const AnalysisDataFrameHeader&               header,
164                                            const std::vector<AnalysisDataValue>&        values,
165                                            const std::vector<AnalysisDataPointSetInfo>& pointSets) :
166     header_(header),
167     values_(values),
168     pointSets_(pointSets)
169 {
170     GMX_ASSERT(!pointSets_.empty(), "There must always be a point set");
171 }
172
173
174 AnalysisDataFrameRef::AnalysisDataFrameRef(const AnalysisDataFrameRef& frame, int firstColumn, int columnCount) :
175     header_(frame.header()),
176     values_(constArrayRefFromArray(&frame.values_[firstColumn], columnCount)),
177     pointSets_(frame.pointSets_)
178 {
179     // FIXME: This doesn't produce a valid internal state, although it does
180     // work in some cases. The point sets cannot be correctly managed here, but
181     // need to be handles by the data proxy class.
182     GMX_ASSERT(firstColumn >= 0, "Invalid first column");
183     GMX_ASSERT(columnCount >= 0, "Invalid column count");
184     GMX_ASSERT(pointSets_.size() == 1U, "Subsets of frames only supported with simple data");
185     GMX_ASSERT(firstColumn + columnCount <= ssize(values_), "Invalid last column");
186 }
187
188
189 bool AnalysisDataFrameRef::allPresent() const
190 {
191     GMX_ASSERT(isValid(), "Invalid data frame accessed");
192     AnalysisDataValuesRef::const_iterator i;
193     for (i = values_.begin(); i != values_.end(); ++i)
194     {
195         if (!i->isPresent())
196         {
197             return false;
198         }
199     }
200     return true;
201 }
202
203 } // namespace gmx