2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,2012,2013, 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.
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.
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.
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.
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.
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.
37 * Implements gmx::AnalysisDataDisplacementModule.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_analysisdata
42 #include "gromacs/analysisdata/modules/displacement.h"
44 #include "gromacs/legacyheaders/smalloc.h"
45 #include "gromacs/legacyheaders/maths.h"
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"
56 /********************************************************************
57 * AnalysisDataDisplacementModule::Impl
61 * Private implementation class for AnalysisDataDisplacementModule.
63 * \ingroup module_analysisdata
65 class AnalysisDataDisplacementModule::Impl
71 //! Maximum number of particles for which the displacements are calculated.
73 //! Maximum time for which the displacements are needed.
75 //! Number of dimensions per data point.
78 //! true if no frames have been read.
80 //! Stores the time of the first frame.
82 //! Stores the time interval between frames.
84 //! Stores the time of the current frame.
86 //! Stores the index in the store for the current positions.
89 //! Maximum number of positions to store for a particle.
91 //! The total number of positions ever stored (can be larger than \p max_store).
95 //! The most recently calculated displacements.
96 std::vector<AnalysisDataValue> currValues_;
98 //! Histogram module for calculating MSD histograms, or NULL if not set.
99 AnalysisDataBinAverageModule *histm;
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),
110 AnalysisDataDisplacementModule::Impl::~Impl()
115 /********************************************************************
116 * AnalysisDataDisplacementModule
119 AnalysisDataDisplacementModule::AnalysisDataDisplacementModule()
126 AnalysisDataDisplacementModule::~AnalysisDataDisplacementModule()
132 AnalysisDataDisplacementModule::setMaxTime(real tmax)
139 AnalysisDataDisplacementModule::setMSDHistogram(
140 AnalysisDataBinAverageModulePointer histm)
142 GMX_RELEASE_ASSERT(_impl->histm == NULL, "Can only set MSD histogram once");
143 _impl->histm = histm.get();
149 AnalysisDataDisplacementModule::tryGetDataFrameInternal(int /*index*/) const
151 return AnalysisDataFrameRef();
156 AnalysisDataDisplacementModule::requestStorageInternal(int /*nframes*/)
163 AnalysisDataDisplacementModule::flags() const
165 return efAllowMulticolumn;
170 AnalysisDataDisplacementModule::dataStarted(AbstractAnalysisData *data)
172 if (data->columnCount() % _impl->ndim != 0)
174 GMX_THROW(APIError("Data has incorrect number of columns"));
176 _impl->nmax = data->columnCount();
177 snew(_impl->oldval, _impl->nmax);
178 _impl->ci = -_impl->nmax;
180 int ncol = _impl->nmax / _impl->ndim + 1;
181 _impl->currValues_.reserve(ncol);
182 setColumnCount(0, ncol);
187 AnalysisDataDisplacementModule::frameStarted(const AnalysisDataFrameHeader &header)
192 _impl->t0 = header.x();
194 else if (_impl->dt <= 0)
196 _impl->dt = header.x() - _impl->t0;
197 if (_impl->dt < 0 || gmx_within_tol(_impl->dt, 0.0, GMX_REAL_EPS))
199 GMX_THROW(APIError("Identical or decreasing frame times"));
204 if (!gmx_within_tol(header.x() - _impl->t, _impl->dt, GMX_REAL_EPS))
206 GMX_THROW(APIError("Frames not evenly spaced"));
209 _impl->t = header.x();
211 // Allocate memory for all the positions once it is possible.
212 if (_impl->max_store == -1 && !_impl->bFirst)
214 _impl->max_store = _impl->nmax * (int)(_impl->tmax/_impl->dt + 1);
215 srenew(_impl->oldval, _impl->max_store);
218 // Increment the index where current positions are stored.
219 _impl->ci += _impl->nmax;
220 if (_impl->ci >= _impl->max_store)
226 for (int i = 0; i < _impl->nmax; ++i)
228 _impl->p[_impl->ci + i].bPres = false;
232 _impl->bFirst = false;
237 AnalysisDataDisplacementModule::pointsAdded(const AnalysisDataPointSetRef &points)
239 if (points.firstColumn() % _impl->ndim != 0
240 || points.columnCount() % _impl->ndim != 0)
242 GMX_THROW(APIError("Partial data points"));
244 for (int i = 0; i < points.columnCount(); ++i)
246 _impl->oldval[_impl->ci + points.firstColumn() + i] = points.y(i);
252 AnalysisDataDisplacementModule::frameFinished(const AnalysisDataFrameHeader & /*header*/)
254 if (_impl->nstored <= 1)
261 if (_impl->nstored == 2)
265 _impl->histm->init(histogramFromBins(0, _impl->max_store / _impl->nmax,
266 _impl->dt).integerBins());
268 moduleManager().notifyDataStart(this);
270 AnalysisDataFrameHeader header(_impl->nstored - 2, _impl->t, 0);
271 moduleManager().notifyFrameStart(header);
273 for (i = _impl->ci - _impl->nmax, step = 1;
274 step < _impl->nstored && i != _impl->ci;
275 i -= _impl->nmax, ++step)
279 i += _impl->max_store;
281 _impl->currValues_.clear();
282 _impl->currValues_.push_back(AnalysisDataValue(step * _impl->dt));
284 for (int j = 0; j < _impl->nmax; j += _impl->ndim, ++k)
288 for (int d = 0; d < _impl->ndim; ++d)
290 real displ = _impl->oldval[_impl->ci + j + d]
291 - _impl->oldval[i + j + d];
292 dist2 += displ * displ;
294 _impl->currValues_.push_back(AnalysisDataValue(dist2));
296 moduleManager().notifyPointsAdd(AnalysisDataPointSetRef(header, _impl->currValues_));
299 moduleManager().notifyFrameFinish(header);
304 AnalysisDataDisplacementModule::dataFinished()
306 if (_impl->nstored >= 2)
308 moduleManager().notifyDataFinish();