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