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