Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / analysisdata / modules / displacement.cpp
blob310572c3df222846da494376ac408b5928a02997
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,2012,2013,2014,2016,2017,2018, 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.
35 /*! \internal \file
36 * \brief
37 * Implements gmx::AnalysisDataDisplacementModule.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_analysisdata
42 #include "gmxpre.h"
44 #include "displacement.h"
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"
54 namespace gmx
57 /********************************************************************
58 * AnalysisDataDisplacementModule::Impl
61 /*! \internal \brief
62 * Private implementation class for AnalysisDataDisplacementModule.
64 * \ingroup module_analysisdata
66 class AnalysisDataDisplacementModule::Impl
68 public:
69 Impl();
70 ~Impl();
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;
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;
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_;
99 //! Histogram module for calculating MSD histograms, or NULL if not set.
100 AnalysisDataBinAverageModule *histm;
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(nullptr),
107 histm(nullptr)
111 AnalysisDataDisplacementModule::Impl::~Impl()
113 sfree(oldval);
116 /********************************************************************
117 * AnalysisDataDisplacementModule
120 AnalysisDataDisplacementModule::AnalysisDataDisplacementModule()
121 : _impl(new Impl())
123 setMultipoint(true);
127 AnalysisDataDisplacementModule::~AnalysisDataDisplacementModule()
132 void
133 AnalysisDataDisplacementModule::setMaxTime(real tmax)
135 _impl->tmax = tmax;
139 void
140 AnalysisDataDisplacementModule::setMSDHistogram(
141 const AnalysisDataBinAverageModulePointer &histm)
143 GMX_RELEASE_ASSERT(_impl->histm == nullptr, "Can only set MSD histogram once");
144 _impl->histm = histm.get();
145 addModule(histm);
149 AnalysisDataFrameRef
150 AnalysisDataDisplacementModule::tryGetDataFrameInternal(int /*index*/) const
152 return AnalysisDataFrameRef();
156 bool
157 AnalysisDataDisplacementModule::requestStorageInternal(int /*nframes*/)
159 return false;
164 AnalysisDataDisplacementModule::flags() const
166 return efAllowMulticolumn;
170 void
171 AnalysisDataDisplacementModule::dataStarted(AbstractAnalysisData *data)
173 if (data->columnCount() % _impl->ndim != 0)
175 GMX_THROW(APIError("Data has incorrect number of columns"));
177 _impl->nmax = data->columnCount();
178 snew(_impl->oldval, _impl->nmax);
179 _impl->ci = -_impl->nmax;
181 int ncol = _impl->nmax / _impl->ndim + 1;
182 _impl->currValues_.reserve(ncol);
183 setColumnCount(0, ncol);
187 void
188 AnalysisDataDisplacementModule::frameStarted(const AnalysisDataFrameHeader &header)
190 // Initialize times.
191 if (_impl->bFirst)
193 _impl->t0 = header.x();
195 else if (_impl->dt <= 0)
197 _impl->dt = header.x() - _impl->t0;
198 if (_impl->dt < 0 || gmx_within_tol(_impl->dt, 0.0, GMX_REAL_EPS))
200 GMX_THROW(APIError("Identical or decreasing frame times"));
203 else
205 if (!gmx_within_tol(header.x() - _impl->t, _impl->dt, GMX_REAL_EPS))
207 GMX_THROW(APIError("Frames not evenly spaced"));
210 _impl->t = header.x();
212 // Allocate memory for all the positions once it is possible.
213 if (_impl->max_store == -1 && !_impl->bFirst)
215 _impl->max_store = _impl->nmax * static_cast<int>(_impl->tmax/_impl->dt + 1);
216 srenew(_impl->oldval, _impl->max_store);
219 // Increment the index where current positions are stored.
220 _impl->ci += _impl->nmax;
221 if (_impl->ci >= _impl->max_store)
223 _impl->ci = 0;
227 for (int i = 0; i < _impl->nmax; ++i)
229 _impl->p[_impl->ci + i].bPres = false;
232 _impl->nstored++;
233 _impl->bFirst = false;
237 void
238 AnalysisDataDisplacementModule::pointsAdded(const AnalysisDataPointSetRef &points)
240 if (points.firstColumn() % _impl->ndim != 0
241 || points.columnCount() % _impl->ndim != 0)
243 GMX_THROW(APIError("Partial data points"));
245 for (int i = 0; i < points.columnCount(); ++i)
247 _impl->oldval[_impl->ci + points.firstColumn() + i] = points.y(i);
252 void
253 AnalysisDataDisplacementModule::frameFinished(const AnalysisDataFrameHeader & /*header*/)
255 if (_impl->nstored <= 1)
257 return;
260 int step, i;
262 if (_impl->nstored == 2)
264 if (_impl->histm)
266 _impl->histm->init(histogramFromBins(0, _impl->max_store / _impl->nmax,
267 _impl->dt).integerBins());
269 moduleManager().notifyDataStart(this);
271 AnalysisDataFrameHeader header(_impl->nstored - 2, _impl->t, 0);
272 moduleManager().notifyFrameStart(header);
274 for (i = _impl->ci - _impl->nmax, step = 1;
275 step < _impl->nstored && i != _impl->ci;
276 i -= _impl->nmax, ++step)
278 if (i < 0)
280 i += _impl->max_store;
282 _impl->currValues_.clear();
283 _impl->currValues_.emplace_back(step * _impl->dt);
284 int k = 1;
285 for (int j = 0; j < _impl->nmax; j += _impl->ndim, ++k)
287 real dist2 = 0.0;
289 for (int d = 0; d < _impl->ndim; ++d)
291 real displ = _impl->oldval[_impl->ci + j + d]
292 - _impl->oldval[i + j + d];
293 dist2 += displ * displ;
295 _impl->currValues_.emplace_back(dist2);
297 moduleManager().notifyPointsAdd(AnalysisDataPointSetRef(header, _impl->currValues_));
300 moduleManager().notifyFrameFinish(header);
304 void
305 AnalysisDataDisplacementModule::dataFinished()
307 if (_impl->nstored >= 2)
309 moduleManager().notifyDataFinish();
313 } // namespace gmx