Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / analysisdata / arraydata.cpp
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11  * Copyright (c) 2001-2009, The GROMACS development team,
12  * check out http://www.gromacs.org for more information.
13
14  * This program is free software; you can redistribute it and/or
15  * modify it under the terms of the GNU General Public License
16  * as published by the Free Software Foundation; either version 2
17  * of the License, or (at your option) any later version.
18  *
19  * If you want to redistribute modifications, please consider that
20  * scientific software is very special. Version control is crucial -
21  * bugs must be traceable. We will be happy to consider code for
22  * inclusion in the official distribution, but derived work must not
23  * be called official GROMACS. Details are found in the README & COPYING
24  * files - if they are missing, get the official version at www.gromacs.org.
25  *
26  * To help us fund GROMACS development, we humbly ask that you cite
27  * the papers on the package - you can find them in the top README file.
28  *
29  * For more info, check our website at http://www.gromacs.org
30  */
31 /*! \internal \file
32  * \brief
33  * Implements classes in arraydata.h.
34  *
35  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
36  * \ingroup module_analysisdata
37  */
38 #include "gromacs/analysisdata/arraydata.h"
39
40 #include <algorithm>
41
42 #include "gromacs/analysisdata/dataframe.h"
43 #include "gromacs/utility/exceptions.h"
44 #include "gromacs/utility/gmxassert.h"
45
46 namespace gmx
47 {
48
49 AbstractAnalysisArrayData::AbstractAnalysisArrayData()
50     : rowCount_(0), xstart_(0.0), xstep_(1.0), bReady_(false)
51 {
52 }
53
54 AbstractAnalysisArrayData::~AbstractAnalysisArrayData()
55 {
56 }
57
58
59 AnalysisDataFrameRef
60 AbstractAnalysisArrayData::tryGetDataFrameInternal(int index) const
61 {
62     if (!isAllocated())
63     {
64         return AnalysisDataFrameRef();
65     }
66     std::vector<AnalysisDataValue>::const_iterator begin
67         = value_.begin() + index * columnCount();
68     return AnalysisDataFrameRef(
69                 AnalysisDataFrameHeader(index, xvalue(index), 0.0),
70                 AnalysisDataValuesRef(begin, begin + columnCount()));
71 }
72
73
74 bool
75 AbstractAnalysisArrayData::requestStorageInternal(int /*nframes*/)
76 {
77     return true;
78 }
79
80
81 void
82 AbstractAnalysisArrayData::setColumnCount(int ncols)
83 {
84     GMX_RELEASE_ASSERT(!isAllocated(),
85                        "Cannot change column count after data has been allocated");
86     AbstractAnalysisData::setColumnCount(ncols);
87 }
88
89
90 void
91 AbstractAnalysisArrayData::setRowCount(int rowCount)
92 {
93     GMX_RELEASE_ASSERT(rowCount > 0, "Invalid number of rows");
94     GMX_RELEASE_ASSERT(!isAllocated(),
95                        "Cannot change row count after data has been allocated");
96     rowCount_ = rowCount;
97 }
98
99
100 void
101 AbstractAnalysisArrayData::allocateValues()
102 {
103     GMX_RELEASE_ASSERT(!isAllocated(), "Can only allocate values once");
104     GMX_RELEASE_ASSERT(rowCount() > 0 && columnCount() > 0,
105                        "Row and column counts must be set before allocating values");
106     value_.resize(rowCount() * columnCount());
107     std::vector<AnalysisDataValue>::iterator i;
108     for (i = value_.begin(); i != value_.end(); ++i)
109     {
110         i->setValue(0.0);
111     }
112 }
113
114
115 void
116 AbstractAnalysisArrayData::setXAxis(real start, real step)
117 {
118     GMX_RELEASE_ASSERT(!bReady_, "X axis cannot be set after data is finished");
119     xstart_ = start;
120     xstep_ = step;
121 }
122
123
124 void
125 AbstractAnalysisArrayData::valuesReady()
126 {
127     GMX_RELEASE_ASSERT(isAllocated(), "There must be some data");
128     if (bReady_)
129     {
130         return;
131     }
132     bReady_ = true;
133
134     std::vector<AnalysisDataValue>::const_iterator valueIter = value_.begin();
135     notifyDataStart();
136     for (int i = 0; i < rowCount(); ++i, valueIter += columnCount())
137     {
138         AnalysisDataFrameHeader header(i, xvalue(i), 0);
139         notifyFrameStart(header);
140         notifyPointsAdd(AnalysisDataPointSetRef(header, 0,
141                             AnalysisDataValuesRef(valueIter,
142                                                   valueIter + columnCount())));
143         notifyFrameFinish(header);
144     }
145     notifyDataFinish();
146 }
147
148
149 void
150 AbstractAnalysisArrayData::copyContents(const AbstractAnalysisArrayData *src,
151                                         AbstractAnalysisArrayData *dest)
152 {
153     GMX_RELEASE_ASSERT(src->isAllocated(), "Source data must not be empty");
154     GMX_RELEASE_ASSERT(!dest->isAllocated(),
155                        "Destination data must not be allocated");
156     dest->setColumnCount(src->columnCount());
157     dest->setRowCount(src->rowCount());
158     dest->allocateValues();
159     dest->setXAxis(src->xstart(), src->xstep());
160     std::copy(src->value_.begin(), src->value_.end(), dest->value_.begin());
161 }
162
163 } // namespace gmx