Merge branch release-4-6
[alexxy/gromacs.git] / src / gromacs / analysisdata / modules / displacement.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2012,2013,2014, 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 gmx::AnalysisDataDisplacementModule.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_analysisdata
41  */
42 #include "gromacs/analysisdata/modules/displacement.h"
43
44 #include "gromacs/legacyheaders/smalloc.h"
45 #include "gromacs/math/utilities.h"
46
47 #include "gromacs/analysisdata/dataframe.h"
48 #include "gromacs/analysisdata/datamodulemanager.h"
49 #include "gromacs/analysisdata/modules/histogram.h"
50 #include "gromacs/utility/exceptions.h"
51 #include "gromacs/utility/gmxassert.h"
52
53 namespace gmx
54 {
55
56 /********************************************************************
57  * AnalysisDataDisplacementModule::Impl
58  */
59
60 /*! \internal \brief
61  * Private implementation class for AnalysisDataDisplacementModule.
62  *
63  * \ingroup module_analysisdata
64  */
65 class AnalysisDataDisplacementModule::Impl
66 {
67     public:
68         Impl();
69         ~Impl();
70
71         //! Maximum number of particles for which the displacements are calculated.
72         int                     nmax;
73         //! Maximum time for which the displacements are needed.
74         real                    tmax;
75         //! Number of dimensions per data point.
76         int                     ndim;
77
78         //! true if no frames have been read.
79         bool                    bFirst;
80         //! Stores the time of the first frame.
81         real                    t0;
82         //! Stores the time interval between frames.
83         real                    dt;
84         //! Stores the time of the current frame.
85         real                    t;
86         //! Stores the index in the store for the current positions.
87         int                     ci;
88
89         //! Maximum number of positions to store for a particle.
90         int                            max_store;
91         //! The total number of positions ever stored (can be larger than \p max_store).
92         int                            nstored;
93         //! Old values.
94         real                          *oldval;
95         //! The most recently calculated displacements.
96         std::vector<AnalysisDataValue> currValues_;
97
98         //! Histogram module for calculating MSD histograms, or NULL if not set.
99         AnalysisDataBinAverageModule *histm;
100 };
101
102 AnalysisDataDisplacementModule::Impl::Impl()
103     : nmax(0), tmax(0.0), ndim(3),
104       bFirst(true), t0(0.0), dt(0.0), t(0.0), ci(0),
105       max_store(-1), nstored(0), oldval(NULL),
106       histm(NULL)
107 {
108 }
109
110 AnalysisDataDisplacementModule::Impl::~Impl()
111 {
112     sfree(oldval);
113 }
114
115 /********************************************************************
116  * AnalysisDataDisplacementModule
117  */
118
119 AnalysisDataDisplacementModule::AnalysisDataDisplacementModule()
120     : _impl(new Impl())
121 {
122     setMultipoint(true);
123 }
124
125
126 AnalysisDataDisplacementModule::~AnalysisDataDisplacementModule()
127 {
128 }
129
130
131 void
132 AnalysisDataDisplacementModule::setMaxTime(real tmax)
133 {
134     _impl->tmax = tmax;
135 }
136
137
138 void
139 AnalysisDataDisplacementModule::setMSDHistogram(
140         AnalysisDataBinAverageModulePointer histm)
141 {
142     GMX_RELEASE_ASSERT(_impl->histm == NULL, "Can only set MSD histogram once");
143     _impl->histm = histm.get();
144     addModule(histm);
145 }
146
147
148 AnalysisDataFrameRef
149 AnalysisDataDisplacementModule::tryGetDataFrameInternal(int /*index*/) const
150 {
151     return AnalysisDataFrameRef();
152 }
153
154
155 bool
156 AnalysisDataDisplacementModule::requestStorageInternal(int /*nframes*/)
157 {
158     return false;
159 }
160
161
162 int
163 AnalysisDataDisplacementModule::flags() const
164 {
165     return efAllowMulticolumn;
166 }
167
168
169 void
170 AnalysisDataDisplacementModule::dataStarted(AbstractAnalysisData *data)
171 {
172     if (data->columnCount() % _impl->ndim != 0)
173     {
174         GMX_THROW(APIError("Data has incorrect number of columns"));
175     }
176     _impl->nmax = data->columnCount();
177     snew(_impl->oldval, _impl->nmax);
178     _impl->ci = -_impl->nmax;
179
180     int ncol = _impl->nmax / _impl->ndim + 1;
181     _impl->currValues_.reserve(ncol);
182     setColumnCount(0, ncol);
183 }
184
185
186 void
187 AnalysisDataDisplacementModule::frameStarted(const AnalysisDataFrameHeader &header)
188 {
189     // Initialize times.
190     if (_impl->bFirst)
191     {
192         _impl->t0 = header.x();
193     }
194     else if (_impl->dt <= 0)
195     {
196         _impl->dt = header.x() - _impl->t0;
197         if (_impl->dt < 0 || gmx_within_tol(_impl->dt, 0.0, GMX_REAL_EPS))
198         {
199             GMX_THROW(APIError("Identical or decreasing frame times"));
200         }
201     }
202     else
203     {
204         if (!gmx_within_tol(header.x() - _impl->t, _impl->dt, GMX_REAL_EPS))
205         {
206             GMX_THROW(APIError("Frames not evenly spaced"));
207         }
208     }
209     _impl->t = header.x();
210
211     // Allocate memory for all the positions once it is possible.
212     if (_impl->max_store == -1 && !_impl->bFirst)
213     {
214         _impl->max_store = _impl->nmax * (int)(_impl->tmax/_impl->dt + 1);
215         srenew(_impl->oldval, _impl->max_store);
216     }
217
218     // Increment the index where current positions are stored.
219     _impl->ci += _impl->nmax;
220     if (_impl->ci >= _impl->max_store)
221     {
222         _impl->ci = 0;
223     }
224
225 /*
226     for (int i = 0; i < _impl->nmax; ++i)
227     {
228         _impl->p[_impl->ci + i].bPres = false;
229     }
230  */
231     _impl->nstored++;
232     _impl->bFirst = false;
233 }
234
235
236 void
237 AnalysisDataDisplacementModule::pointsAdded(const AnalysisDataPointSetRef &points)
238 {
239     if (points.firstColumn() % _impl->ndim != 0
240         || points.columnCount() % _impl->ndim != 0)
241     {
242         GMX_THROW(APIError("Partial data points"));
243     }
244     for (int i = 0; i < points.columnCount(); ++i)
245     {
246         _impl->oldval[_impl->ci + points.firstColumn() + i] = points.y(i);
247     }
248 }
249
250
251 void
252 AnalysisDataDisplacementModule::frameFinished(const AnalysisDataFrameHeader & /*header*/)
253 {
254     if (_impl->nstored <= 1)
255     {
256         return;
257     }
258
259     int step, i;
260
261     if (_impl->nstored == 2)
262     {
263         if (_impl->histm)
264         {
265             _impl->histm->init(histogramFromBins(0, _impl->max_store / _impl->nmax,
266                                                  _impl->dt).integerBins());
267         }
268         moduleManager().notifyDataStart(this);
269     }
270     AnalysisDataFrameHeader header(_impl->nstored - 2, _impl->t, 0);
271     moduleManager().notifyFrameStart(header);
272
273     for (i = _impl->ci - _impl->nmax, step = 1;
274          step < _impl->nstored && i != _impl->ci;
275          i -= _impl->nmax, ++step)
276     {
277         if (i < 0)
278         {
279             i += _impl->max_store;
280         }
281         _impl->currValues_.clear();
282         _impl->currValues_.push_back(AnalysisDataValue(step * _impl->dt));
283         int k = 1;
284         for (int j = 0; j < _impl->nmax; j += _impl->ndim, ++k)
285         {
286             real dist2 = 0.0;
287
288             for (int d = 0; d < _impl->ndim; ++d)
289             {
290                 real displ = _impl->oldval[_impl->ci + j + d]
291                     - _impl->oldval[i + j + d];
292                 dist2 += displ * displ;
293             }
294             _impl->currValues_.push_back(AnalysisDataValue(dist2));
295         }
296         moduleManager().notifyPointsAdd(AnalysisDataPointSetRef(header, _impl->currValues_));
297     }
298
299     moduleManager().notifyFrameFinish(header);
300 }
301
302
303 void
304 AnalysisDataDisplacementModule::dataFinished()
305 {
306     if (_impl->nstored >= 2)
307     {
308         moduleManager().notifyDataFinish();
309     }
310 }
311
312 } // namespace gmx