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