2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,2012,2013,2014,2015,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.
37 * Declares analysis data modules for calculating histograms.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
41 * \ingroup module_analysisdata
43 #ifndef GMX_ANALYSISDATA_MODULES_HISTOGRAM_H
44 #define GMX_ANALYSISDATA_MODULES_HISTOGRAM_H
48 #include "gromacs/analysisdata/abstractdata.h"
49 #include "gromacs/analysisdata/arraydata.h"
50 #include "gromacs/analysisdata/datamodule.h"
55 class AnalysisHistogramSettings
;
58 * Provides "named parameter" idiom for constructing histograms.
60 * \see histogramFromBins()
61 * \see histogramFromRange()
63 * Methods in this class do not throw.
66 * \ingroup module_analysisdata
68 class AnalysisHistogramSettingsInitializer
72 * Creates an empty initializer.
74 * Should not be called directly, but histogramFromRange() or
75 * histogramFromBins() should be used instead.
77 AnalysisHistogramSettingsInitializer();
80 * Sets the first bin location.
82 * Typically should not be called directly, but through
83 * histogramFromBins().
85 AnalysisHistogramSettingsInitializer
&start(real min
)
86 { min_
= min
; return *this; }
88 * Sets the number of bins in the histogram.
90 * If only the first bin location is specified, this value is required
91 * (and automatically provided if histogramFromBins() is used).
92 * If both the first and last bins are specified, either this value or
93 * binWidth() is required.
95 AnalysisHistogramSettingsInitializer
&binCount(int binCount
)
96 { binCount_
= binCount
; return *this; }
98 * Sets the first and last bin locations.
100 * Typically should not be called directly, but through
101 * histogramFromRange().
103 AnalysisHistogramSettingsInitializer
&range(real min
, real max
)
104 { min_
= min
; max_
= max
; return *this; }
106 * Sets the bin width of the histogram.
108 * If only the first bin location is specified, this value is required
109 * (and automatically provided if histogramFromBins() is used).
110 * If both the first and last bins are specified, either this value or
111 * binCount() is required.
112 * If a bin width is provided with both first and last bin locations,
113 * and the given bin width does not divide the range exactly, the last
114 * bin location is adjusted to match.
116 AnalysisHistogramSettingsInitializer
&binWidth(real binWidth
)
117 { binWidth_
= binWidth
; return *this; }
119 * Indicate that first and last bin locations to specify bin centers.
121 * If set, the first and last bin locations are interpreted as bin
123 * If not set (the default), the first and last bin locations are
124 * interpreted as the edges of the whole histogram.
126 * Cannot be specified together with roundRange().
128 AnalysisHistogramSettingsInitializer
&integerBins(bool enabled
= true)
129 { bIntegerBins_
= enabled
; return *this; }
131 * Round first and last bin locations.
133 * If set, the resulting histogram will cover the range specified, but
134 * the actual bin locations will be rounded such that the edges fall
135 * on multiples of the bin width.
136 * Only implemented when both first and last bin location and bin width
138 * Cannot be specified together with integerBins() or with binCount().
140 AnalysisHistogramSettingsInitializer
&roundRange(bool enabled
= true)
141 { bRoundRange_
= enabled
; return *this; }
143 * Sets the histogram to match all values.
145 * If set, the histogram behaves as if the bins at the ends extended to
148 AnalysisHistogramSettingsInitializer
&includeAll(bool enabled
= true)
149 { bIncludeAll_
= enabled
; return *this; }
160 friend class AnalysisHistogramSettings
;
164 * Initializes a histogram using a range and a bin width.
170 inline AnalysisHistogramSettingsInitializer
171 histogramFromRange(real min
, real max
)
173 return AnalysisHistogramSettingsInitializer().range(min
, max
);
177 * Initializes a histogram using bin width and the number of bins.
183 inline AnalysisHistogramSettingsInitializer
184 histogramFromBins(real start
, int nbins
, real binwidth
)
186 return AnalysisHistogramSettingsInitializer()
187 .start(start
).binCount(nbins
).binWidth(binwidth
);
192 * Contains parameters that specify histogram bin locations.
194 * Methods in this class do not throw.
197 * \ingroup module_analysisdata
199 class AnalysisHistogramSettings
202 //! Initializes undefined parameters.
203 AnalysisHistogramSettings();
205 * Initializes parameters based on a named parameter object.
207 * This constructor is not explicit to allow initialization of
208 * histograms directly from AnalysisHistogramSettingsInitializer:
210 gmx::AnalysisDataSimpleHistogramModule *hist =
211 new gmx::AnalysisDataSimpleHistogramModule(
212 histogramFromRange(0.0, 5.0).binWidth(0.5));
215 AnalysisHistogramSettings(const AnalysisHistogramSettingsInitializer
&settings
);
217 //! Returns the left edge of the first bin.
218 real
firstEdge() const { return firstEdge_
; }
219 //! Returns the right edge of the first bin.
220 real
lastEdge() const { return lastEdge_
; }
221 //! Returns the number of bins in the histogram.
222 int binCount() const { return binCount_
; }
223 //! Returns the width of a bin in the histogram.
224 real
binWidth() const { return binWidth_
; }
225 //! Whether values beyond the edges are mapped to the edge bins.
226 bool includeAll() const { return bAll_
; }
227 //! Returns a zero-based bin index for a value, or -1 if not in range.
228 int findBin(real y
) const;
234 real inverseBinWidth_
;
240 class AbstractAverageHistogram
;
242 //! Smart pointer to manage an AbstractAverageHistogram object.
243 typedef std::unique_ptr
<AbstractAverageHistogram
> AverageHistogramPointer
;
246 * Base class for representing histograms averaged over frames.
248 * The averaging module for a per-frame histogram is always created by the
249 * histogram module class (e.g., AnalysisDataSimpleHistogramModule), and can be
250 * accessed using, e.g., AnalysisDataSimpleHistogramModule::averager().
251 * The user can alter some properties of the average histogram directly, but
252 * the main use of the object is to postprocess the histogram once the
253 * calculation is finished.
255 * This class can represent multiple histograms in one object: each column in
256 * the data is an independent histogram.
257 * The X values correspond to center of the bins, except for a cumulative
258 * histogram made with makeCumulative().
261 * \ingroup module_analysisdata
263 class AbstractAverageHistogram
: public AbstractAnalysisArrayData
266 ~AbstractAverageHistogram() override
;
268 //! Returns bin properties for the histogram.
269 const AnalysisHistogramSettings
&settings() const { return settings_
; }
272 * Creates a copy of the histogram with double the bin width.
274 * \param[in] bIntegerBins If `true`, the first bin in the result will
275 * cover the first bin from the source. Otherwise, the first bin
276 * will cover first two bins from the source.
277 * \throws std::bad_alloc if out of memory.
279 * The caller is responsible of deleting the returned object.
281 AverageHistogramPointer
resampleDoubleBinWidth(bool bIntegerBins
) const;
283 * Creates a deep copy of the histogram.
285 * \throws std::bad_alloc if out of memory.
287 * The returned histogram is not necessarily of the same dynamic type
288 * as the original object, but contains the same data from the point of
289 * view of the AbstractAverageHistogram interface.
291 * The caller is responsible of deleting the returned object.
293 AverageHistogramPointer
clone() const;
294 //! Normalizes the histogram such that the integral over it is one.
295 void normalizeProbability();
297 * Makes the histograms cumulative by summing up each bin to all bins
300 * The X values in the data are adjusted such that they match the right
301 * edges of bins instead of bin centers.
303 void makeCumulative();
304 //! Scales a single histogram by a uniform scaling factor.
305 void scaleSingle(int index
, real factor
);
306 //! Scales all histograms by a uniform scaling factor.
307 void scaleAll(real factor
);
308 //! Scales the value of each bin by a different scaling factor.
309 void scaleAllByVector(const real factor
[]);
311 * Notifies attached modules of the histogram data.
313 * After this function has been called, it is no longer possible to
314 * alter the histogram.
316 void done() { AbstractAnalysisArrayData::valuesReady(); }
320 * Creates a histogram module with undefined bins.
322 * Bin parameters must be defined with init() before data input is
325 AbstractAverageHistogram();
326 //! Creates a histogram module with defined bin parameters.
327 explicit AbstractAverageHistogram(const AnalysisHistogramSettings
&settings
);
330 * (Re)initializes the histogram from settings.
332 void init(const AnalysisHistogramSettings
&settings
);
335 AnalysisHistogramSettings settings_
;
337 // Copy and assign disallowed by base.
342 * Data module for per-frame histograms.
344 * Output data contains the same number of frames and data sets as the input
345 * data. Each frame contains the histogram(s) for the points in that frame.
346 * Each input data set is processed independently into the corresponding output
347 * data set. Missing values are ignored.
348 * All input columns for a data set are averaged into the same histogram.
349 * The number of columns for all data sets equals the number of bins in the
352 * The histograms are accumulated as 64-bit integers within a frame and summed
353 * in double precision across frames, even if the output data is in single
357 * \ingroup module_analysisdata
359 class AnalysisDataSimpleHistogramModule
: public AbstractAnalysisData
,
360 public AnalysisDataModuleParallel
364 * Creates a histogram module with undefined bins.
366 * Bin parameters must be defined with init() before data input is
369 AnalysisDataSimpleHistogramModule();
370 //! Creates a histogram module with defined bin parameters.
371 explicit AnalysisDataSimpleHistogramModule(const AnalysisHistogramSettings
&settings
);
372 ~AnalysisDataSimpleHistogramModule() override
;
375 * (Re)initializes the histogram from settings.
377 void init(const AnalysisHistogramSettings
&settings
);
380 * Returns the average histogram over all frames.
382 * Can be called already before the histogram is calculated to
383 * customize the way the average histogram is calculated.
385 * \see AbstractAverageHistogram
387 AbstractAverageHistogram
&averager();
389 //! Returns bin properties for the histogram.
390 const AnalysisHistogramSettings
&settings() const;
392 int frameCount() const override
;
394 int flags() const override
;
396 bool parallelDataStarted(
397 AbstractAnalysisData
*data
,
398 const AnalysisDataParallelOptions
&options
) override
;
399 void frameStarted(const AnalysisDataFrameHeader
&header
) override
;
400 void pointsAdded(const AnalysisDataPointSetRef
&points
) override
;
401 void frameFinished(const AnalysisDataFrameHeader
&header
) override
;
402 void frameFinishedSerial(int frameIndex
) override
;
403 void dataFinished() override
;
406 AnalysisDataFrameRef
tryGetDataFrameInternal(int index
) const override
;
407 bool requestStorageInternal(int nframes
) override
;
411 PrivateImplPointer
<Impl
> impl_
;
413 // Copy and assign disallowed by base.
418 * Data module for per-frame weighted histograms.
420 * Output data contains the same number of frames and data sets as the input
421 * data. Each frame contains the histogram(s) for the points in that frame,
422 * interpreted such that the first column passed to pointsAdded() determines
423 * the bin and the rest give weights to be added to that bin (input data should
424 * have at least two colums, and at least two columns should be added at the
426 * Each input data set is processed independently into the corresponding output
428 * All input columns for a data set are averaged into the same histogram.
429 * The number of columns for all data sets equals the number of bins in the
432 * The histograms are accumulated in double precision, even if the output data
433 * is in single precision.
436 * \ingroup module_analysisdata
438 class AnalysisDataWeightedHistogramModule
: public AbstractAnalysisData
,
439 public AnalysisDataModuleParallel
442 //! \copydoc AnalysisDataSimpleHistogramModule::AnalysisDataSimpleHistogramModule()
443 AnalysisDataWeightedHistogramModule();
444 //! \copydoc AnalysisDataSimpleHistogramModule::AnalysisDataSimpleHistogramModule(const AnalysisHistogramSettings &)
445 explicit AnalysisDataWeightedHistogramModule(const AnalysisHistogramSettings
&settings
);
446 ~AnalysisDataWeightedHistogramModule() override
;
448 //! \copydoc AnalysisDataSimpleHistogramModule::init()
449 void init(const AnalysisHistogramSettings
&settings
);
451 //! \copydoc AnalysisDataSimpleHistogramModule::averager()
452 AbstractAverageHistogram
&averager();
454 //! \copydoc AnalysisDataSimpleHistogramModule::settings()
455 const AnalysisHistogramSettings
&settings() const;
457 int frameCount() const override
;
459 int flags() const override
;
461 bool parallelDataStarted(
462 AbstractAnalysisData
*data
,
463 const AnalysisDataParallelOptions
&options
) override
;
464 void frameStarted(const AnalysisDataFrameHeader
&header
) override
;
465 void pointsAdded(const AnalysisDataPointSetRef
&points
) override
;
466 void frameFinished(const AnalysisDataFrameHeader
&header
) override
;
467 void frameFinishedSerial(int frameIndex
) override
;
468 void dataFinished() override
;
471 AnalysisDataFrameRef
tryGetDataFrameInternal(int index
) const override
;
472 bool requestStorageInternal(int nframes
) override
;
476 PrivateImplPointer
<Impl
> impl_
;
478 // Copy and assign disallowed by base.
483 * Data module for bin averages.
485 * Output data contains one row for each bin; see AbstractAverageHistogram.
486 * Output data contains one column for each input data set.
487 * The value in a column is the average over all frames of that data set for
489 * The input data is interpreted such that the first column passed to
490 * pointsAdded() determines the bin and the rest give values to be added to
491 * that bin (input data should have at least two colums, and at least two
492 * columns should be added at the same time).
493 * All input columns for a data set are averaged into the same histogram.
496 * \ingroup module_analysisdata
498 class AnalysisDataBinAverageModule
: public AbstractAnalysisArrayData
,
499 public AnalysisDataModuleSerial
502 //! \copydoc AnalysisDataSimpleHistogramModule::AnalysisDataSimpleHistogramModule()
503 AnalysisDataBinAverageModule();
504 //! \copydoc AnalysisDataSimpleHistogramModule::AnalysisDataSimpleHistogramModule(const AnalysisHistogramSettings &)
505 explicit AnalysisDataBinAverageModule(const AnalysisHistogramSettings
&settings
);
506 ~AnalysisDataBinAverageModule() override
;
508 //! \copydoc AnalysisDataSimpleHistogramModule::init()
509 void init(const AnalysisHistogramSettings
&settings
);
511 //! \copydoc AnalysisDataSimpleHistogramModule::settings()
512 const AnalysisHistogramSettings
&settings() const;
514 int flags() const override
;
516 void dataStarted(AbstractAnalysisData
*data
) override
;
517 void frameStarted(const AnalysisDataFrameHeader
&header
) override
;
518 void pointsAdded(const AnalysisDataPointSetRef
&points
) override
;
519 void frameFinished(const AnalysisDataFrameHeader
&header
) override
;
520 void dataFinished() override
;
525 PrivateImplPointer
<Impl
> impl_
;
527 // Copy and assign disallowed by base.
530 //! Smart pointer to manage an AnalysisDataSimpleHistogramModule object.
531 typedef std::shared_ptr
<AnalysisDataSimpleHistogramModule
>
532 AnalysisDataSimpleHistogramModulePointer
;
533 //! Smart pointer to manage an AnalysisDataWeightedHistogramModule object.
534 typedef std::shared_ptr
<AnalysisDataWeightedHistogramModule
>
535 AnalysisDataWeightedHistogramModulePointer
;
536 //! Smart pointer to manage an AnalysisDataBinAverageModule object.
537 typedef std::shared_ptr
<AnalysisDataBinAverageModule
>
538 AnalysisDataBinAverageModulePointer
;