bea6e7ad33042fcd6e305f9547ffd9436a57cfe7
[alexxy/gromacs.git] / src / gromacs / analysisdata / arraydata.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2012,2013,2014,2015,2017, 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 arraydata.h.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_analysisdata
41  */
42 #include "gmxpre.h"
43
44 #include "arraydata.h"
45
46 #include <algorithm>
47
48 #include "gromacs/analysisdata/dataframe.h"
49 #include "gromacs/analysisdata/datamodulemanager.h"
50 #include "gromacs/utility/exceptions.h"
51 #include "gromacs/utility/gmxassert.h"
52
53 namespace gmx
54 {
55
56 AbstractAnalysisArrayData::AbstractAnalysisArrayData()
57     : rowCount_(0), pointSetInfo_(0, 0, 0, 0), xstep_(1.0),
58       bUniformX_(true), bReady_(false)
59 {
60     xvalue_.push_back(0);
61 }
62
63 AbstractAnalysisArrayData::~AbstractAnalysisArrayData()
64 {
65 }
66
67
68 AnalysisDataFrameRef
69 AbstractAnalysisArrayData::tryGetDataFrameInternal(int index) const
70 {
71     if (!isAllocated())
72     {
73         return AnalysisDataFrameRef();
74     }
75     return AnalysisDataFrameRef(
76             AnalysisDataFrameHeader(index, xvalue(index), 0.0),
77             makeConstArrayRef(value_).
78                 subArray(index * columnCount(), columnCount()),
79             constArrayRefFromArray(&pointSetInfo_, 1));
80 }
81
82
83 bool
84 AbstractAnalysisArrayData::requestStorageInternal(int /*nframes*/)
85 {
86     return true;
87 }
88
89
90 void
91 AbstractAnalysisArrayData::setColumnCount(int ncols)
92 {
93     GMX_RELEASE_ASSERT(!isAllocated(),
94                        "Cannot change column count after data has been allocated");
95     AbstractAnalysisData::setColumnCount(0, ncols);
96     pointSetInfo_ = AnalysisDataPointSetInfo(0, ncols, 0, 0);
97 }
98
99
100 void
101 AbstractAnalysisArrayData::setRowCount(int rowCount)
102 {
103     GMX_RELEASE_ASSERT(rowCount > 0, "Invalid number of rows");
104     GMX_RELEASE_ASSERT(!isAllocated(),
105                        "Cannot change row count after data has been allocated");
106     GMX_RELEASE_ASSERT(bUniformX_ || xvalue_.empty()
107                        || rowCount == static_cast<int>(xvalue_.size()),
108                        "X axis set with setXAxisValue() does not match the row count");
109     xvalue_.resize(rowCount);
110     if (bUniformX_ && rowCount > rowCount_)
111     {
112         for (int i = rowCount_; i < rowCount; ++i)
113         {
114             xvalue_[i] = xvalue_[0] + i * xstep_;
115         }
116     }
117     rowCount_ = rowCount;
118 }
119
120
121 void
122 AbstractAnalysisArrayData::allocateValues()
123 {
124     GMX_RELEASE_ASSERT(!isAllocated(), "Can only allocate values once");
125     GMX_RELEASE_ASSERT(rowCount() > 0 && columnCount() > 0,
126                        "Row and column counts must be set before allocating values");
127     value_.resize(rowCount() * columnCount());
128     std::vector<AnalysisDataValue>::iterator i;
129     for (i = value_.begin(); i != value_.end(); ++i)
130     {
131         i->setValue(0.0);
132     }
133 }
134
135
136 void
137 AbstractAnalysisArrayData::setXAxis(real start, real step)
138 {
139     GMX_RELEASE_ASSERT(!bReady_, "X axis cannot be set after data is finished");
140     xvalue_[0] = start;
141     xstep_     = step;
142     bUniformX_ = true;
143     for (int i = 0; i < rowCount_; ++i)
144     {
145         xvalue_[i] = start + i * xstep_;
146     }
147 }
148
149
150 void
151 AbstractAnalysisArrayData::setXAxisValue(int row, real value)
152 {
153     GMX_RELEASE_ASSERT(!bReady_, "X axis cannot be set after data is finished");
154     if (rowCount_ > 0)
155     {
156         GMX_RELEASE_ASSERT(row >= 0 && row < rowCount(), "Row index out of range");
157     }
158     else if (row >= static_cast<int>(xvalue_.size()))
159     {
160         xvalue_.resize(row + 1);
161     }
162     bUniformX_   = false;
163     xstep_       = 0.0;
164     xvalue_[row] = value;
165 }
166
167
168 void
169 AbstractAnalysisArrayData::valuesReady()
170 {
171     GMX_RELEASE_ASSERT(isAllocated(), "There must be some data");
172     if (bReady_)
173     {
174         return;
175     }
176     bReady_ = true;
177
178     AnalysisDataModuleManager                     &modules   = moduleManager();
179     modules.notifyDataStart(this);
180     for (int i = 0; i < rowCount(); ++i)
181     {
182         AnalysisDataFrameHeader header(i, xvalue(i), 0);
183         modules.notifyFrameStart(header);
184         modules.notifyPointsAdd(
185                 AnalysisDataPointSetRef(
186                         header, pointSetInfo_,
187                         makeConstArrayRef(value_).
188                             subArray(i*columnCount(), columnCount())));
189         modules.notifyFrameFinish(header);
190     }
191     modules.notifyDataFinish();
192 }
193
194
195 void
196 AbstractAnalysisArrayData::copyContents(const AbstractAnalysisArrayData *src,
197                                         AbstractAnalysisArrayData       *dest)
198 {
199     GMX_RELEASE_ASSERT(src->isAllocated(), "Source data must not be empty");
200     GMX_RELEASE_ASSERT(!dest->isAllocated(),
201                        "Destination data must not be allocated");
202     dest->setColumnCount(src->columnCount());
203     dest->setRowCount(src->rowCount());
204     dest->allocateValues();
205     dest->xstep_     = src->xstep_;
206     dest->bUniformX_ = src->bUniformX_;
207     std::copy(src->xvalue_.begin(), src->xvalue_.end(), dest->xvalue_.begin());
208     std::copy(src->value_.begin(), src->value_.end(), dest->value_.begin());
209 }
210
211 } // namespace gmx