7a2f8c5ec92fe486969a7cbe9a5949a6f5078cf4
[alexxy/gromacs.git] / src / gromacs / analysisdata / modules / histogram.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2012,2013,2014,2015,2018,2019, 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 /*! \file
36  * \brief
37  * Declares analysis data modules for calculating histograms.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \inpublicapi
41  * \ingroup module_analysisdata
42  */
43 #ifndef GMX_ANALYSISDATA_MODULES_HISTOGRAM_H
44 #define GMX_ANALYSISDATA_MODULES_HISTOGRAM_H
45
46 #include <memory>
47
48 #include "gromacs/analysisdata/abstractdata.h"
49 #include "gromacs/analysisdata/arraydata.h"
50 #include "gromacs/analysisdata/datamodule.h"
51
52 namespace gmx
53 {
54
55 class AnalysisHistogramSettings;
56
57 /*! \brief
58  * Provides "named parameter" idiom for constructing histograms.
59  *
60  * \see histogramFromBins()
61  * \see histogramFromRange()
62  *
63  * Methods in this class do not throw.
64  *
65  * \inpublicapi
66  * \ingroup module_analysisdata
67  */
68 class AnalysisHistogramSettingsInitializer
69 {
70 public:
71     /*! \brief
72      * Creates an empty initializer.
73      *
74      * Should not be called directly, but histogramFromRange() or
75      * histogramFromBins() should be used instead.
76      */
77     AnalysisHistogramSettingsInitializer();
78
79     /*! \brief
80      * Sets the first bin location.
81      *
82      * Typically should not be called directly, but through
83      * histogramFromBins().
84      */
85     AnalysisHistogramSettingsInitializer& start(real min)
86     {
87         min_ = min;
88         return *this;
89     }
90     /*! \brief
91      * Sets the number of bins in the histogram.
92      *
93      * If only the first bin location is specified, this value is required
94      * (and automatically provided if histogramFromBins() is used).
95      * If both the first and last bins are specified, either this value or
96      * binWidth() is required.
97      */
98     AnalysisHistogramSettingsInitializer& binCount(int binCount)
99     {
100         binCount_ = binCount;
101         return *this;
102     }
103     /*! \brief
104      * Sets the first and last bin locations.
105      *
106      * Typically should not be called directly, but through
107      * histogramFromRange().
108      */
109     AnalysisHistogramSettingsInitializer& range(real min, real max)
110     {
111         min_ = min;
112         max_ = max;
113         return *this;
114     }
115     /*! \brief
116      * Sets the bin width of the histogram.
117      *
118      * If only the first bin location is specified, this value is required
119      * (and automatically provided if histogramFromBins() is used).
120      * If both the first and last bins are specified, either this value or
121      * binCount() is required.
122      * If a bin width is provided with both first and last bin locations,
123      * and the given bin width does not divide the range exactly, the last
124      * bin location is adjusted to match.
125      */
126     AnalysisHistogramSettingsInitializer& binWidth(real binWidth)
127     {
128         binWidth_ = binWidth;
129         return *this;
130     }
131     /*! \brief
132      * Indicate that first and last bin locations to specify bin centers.
133      *
134      * If set, the first and last bin locations are interpreted as bin
135      * centers.
136      * If not set (the default), the first and last bin locations are
137      * interpreted as the edges of the whole histogram.
138      *
139      * Cannot be specified together with roundRange().
140      */
141     AnalysisHistogramSettingsInitializer& integerBins(bool enabled = true)
142     {
143         bIntegerBins_ = enabled;
144         return *this;
145     }
146     /*! \brief
147      * Round first and last bin locations.
148      *
149      * If set, the resulting histogram will cover the range specified, but
150      * the actual bin locations will be rounded such that the edges fall
151      * on multiples of the bin width.
152      * Only implemented when both first and last bin location and bin width
153      * are defined.
154      * Cannot be specified together with integerBins() or with binCount().
155      */
156     AnalysisHistogramSettingsInitializer& roundRange(bool enabled = true)
157     {
158         bRoundRange_ = enabled;
159         return *this;
160     }
161     /*! \brief
162      * Sets the histogram to match all values.
163      *
164      * If set, the histogram behaves as if the bins at the ends extended to
165      * +-infinity.
166      */
167     AnalysisHistogramSettingsInitializer& includeAll(bool enabled = true)
168     {
169         bIncludeAll_ = enabled;
170         return *this;
171     }
172
173 private:
174     real min_;
175     real max_;
176     real binWidth_;
177     int  binCount_;
178     bool bIntegerBins_;
179     bool bRoundRange_;
180     bool bIncludeAll_;
181
182     friend class AnalysisHistogramSettings;
183 };
184
185 /*! \brief
186  * Initializes a histogram using a range and a bin width.
187  *
188  * Does not throw.
189  *
190  * \inpublicapi
191  */
192 inline AnalysisHistogramSettingsInitializer histogramFromRange(real min, real max)
193 {
194     return AnalysisHistogramSettingsInitializer().range(min, max);
195 }
196
197 /*! \brief
198  * Initializes a histogram using bin width and the number of bins.
199  *
200  * Does not throw.
201  *
202  * \inpublicapi
203  */
204 inline AnalysisHistogramSettingsInitializer histogramFromBins(real start, int nbins, real binwidth)
205 {
206     return AnalysisHistogramSettingsInitializer().start(start).binCount(nbins).binWidth(binwidth);
207 }
208
209
210 /*! \brief
211  * Contains parameters that specify histogram bin locations.
212  *
213  * Methods in this class do not throw.
214  *
215  * \inpublicapi
216  * \ingroup module_analysisdata
217  */
218 class AnalysisHistogramSettings
219 {
220 public:
221     //! Initializes undefined parameters.
222     AnalysisHistogramSettings();
223     /*! \brief
224      * Initializes parameters based on a named parameter object.
225      *
226      * This constructor is not explicit to allow initialization of
227      * histograms directly from AnalysisHistogramSettingsInitializer:
228      * \code
229        gmx::AnalysisDataSimpleHistogramModule *hist =
230                new gmx::AnalysisDataSimpleHistogramModule(
231                        histogramFromRange(0.0, 5.0).binWidth(0.5));
232      * \endcode
233      */
234     AnalysisHistogramSettings(const AnalysisHistogramSettingsInitializer& settings);
235
236     //! Returns the left edge of the first bin.
237     real firstEdge() const { return firstEdge_; }
238     //! Returns the right edge of the first bin.
239     real lastEdge() const { return lastEdge_; }
240     //! Returns the number of bins in the histogram.
241     int binCount() const { return binCount_; }
242     //! Returns the width of a bin in the histogram.
243     real binWidth() const { return binWidth_; }
244     //! Whether values beyond the edges are mapped to the edge bins.
245     bool includeAll() const { return bAll_; }
246     //! Returns a zero-based bin index for a value, or -1 if not in range.
247     int findBin(real y) const;
248
249 private:
250     real firstEdge_;
251     real lastEdge_;
252     real binWidth_;
253     real inverseBinWidth_;
254     int  binCount_;
255     bool bAll_;
256 };
257
258
259 class AbstractAverageHistogram;
260
261 //! Smart pointer to manage an AbstractAverageHistogram object.
262 typedef std::unique_ptr<AbstractAverageHistogram> AverageHistogramPointer;
263
264 /*! \brief
265  * Base class for representing histograms averaged over frames.
266  *
267  * The averaging module for a per-frame histogram is always created by the
268  * histogram module class (e.g., AnalysisDataSimpleHistogramModule), and can be
269  * accessed using, e.g., AnalysisDataSimpleHistogramModule::averager().
270  * The user can alter some properties of the average histogram directly, but
271  * the main use of the object is to postprocess the histogram once the
272  * calculation is finished.
273  *
274  * This class can represent multiple histograms in one object: each column in
275  * the data is an independent histogram.
276  * The X values correspond to center of the bins, except for a cumulative
277  * histogram made with makeCumulative().
278  *
279  * \inpublicapi
280  * \ingroup module_analysisdata
281  */
282 class AbstractAverageHistogram : public AbstractAnalysisArrayData
283 {
284 public:
285     ~AbstractAverageHistogram() override;
286
287     //! Returns bin properties for the histogram.
288     const AnalysisHistogramSettings& settings() const { return settings_; }
289
290     /*! \brief
291      * Creates a copy of the histogram with double the bin width.
292      *
293      * \param[in] bIntegerBins If `true`, the first bin in the result will
294      *     cover the first bin from the source. Otherwise, the first bin
295      *     will cover first two bins from the source.
296      * \throws std::bad_alloc if out of memory.
297      *
298      * The caller is responsible of deleting the returned object.
299      */
300     AverageHistogramPointer resampleDoubleBinWidth(bool bIntegerBins) const;
301     /*! \brief
302      * Creates a deep copy of the histogram.
303      *
304      * \throws std::bad_alloc if out of memory.
305      *
306      * The returned histogram is not necessarily of the same dynamic type
307      * as the original object, but contains the same data from the point of
308      * view of the AbstractAverageHistogram interface.
309      *
310      * The caller is responsible of deleting the returned object.
311      */
312     AverageHistogramPointer clone() const;
313     //! Normalizes the histogram such that the integral over it is one.
314     void normalizeProbability();
315     /*! \brief
316      * Makes the histograms cumulative by summing up each bin to all bins
317      * after it.
318      *
319      * The X values in the data are adjusted such that they match the right
320      * edges of bins instead of bin centers.
321      */
322     void makeCumulative();
323     //! Scales a single histogram by a uniform scaling factor.
324     void scaleSingle(int index, real factor);
325     //! Scales all histograms by a uniform scaling factor.
326     void scaleAll(real factor);
327     //! Scales the value of each bin by a different scaling factor.
328     void scaleAllByVector(const real factor[]);
329     /*! \brief
330      * Notifies attached modules of the histogram data.
331      *
332      * After this function has been called, it is no longer possible to
333      * alter the histogram.
334      */
335     void done() { AbstractAnalysisArrayData::valuesReady(); }
336
337 protected:
338     /*! \brief
339      * Creates a histogram module with undefined bins.
340      *
341      * Bin parameters must be defined with init() before data input is
342      * started.
343      */
344     AbstractAverageHistogram();
345     //! Creates a histogram module with defined bin parameters.
346     explicit AbstractAverageHistogram(const AnalysisHistogramSettings& settings);
347
348     /*! \brief
349      * (Re)initializes the histogram from settings.
350      */
351     void init(const AnalysisHistogramSettings& settings);
352
353 private:
354     AnalysisHistogramSettings settings_;
355
356     // Copy and assign disallowed by base.
357 };
358
359
360 /*! \brief
361  * Data module for per-frame histograms.
362  *
363  * Output data contains the same number of frames and data sets as the input
364  * data.  Each frame contains the histogram(s) for the points in that frame.
365  * Each input data set is processed independently into the corresponding output
366  * data set.  Missing values are ignored.
367  * All input columns for a data set are averaged into the same histogram.
368  * The number of columns for all data sets equals the number of bins in the
369  * histogram.
370  *
371  * The histograms are accumulated as 64-bit integers within a frame and summed
372  * in double precision across frames, even if the output data is in single
373  * precision.
374  *
375  * \inpublicapi
376  * \ingroup module_analysisdata
377  */
378 class AnalysisDataSimpleHistogramModule : public AbstractAnalysisData, public AnalysisDataModuleParallel
379 {
380 public:
381     /*! \brief
382      * Creates a histogram module with undefined bins.
383      *
384      * Bin parameters must be defined with init() before data input is
385      * started.
386      */
387     AnalysisDataSimpleHistogramModule();
388     //! Creates a histogram module with defined bin parameters.
389     explicit AnalysisDataSimpleHistogramModule(const AnalysisHistogramSettings& settings);
390     ~AnalysisDataSimpleHistogramModule() override;
391
392     /*! \brief
393      * (Re)initializes the histogram from settings.
394      */
395     void init(const AnalysisHistogramSettings& settings);
396
397     /*! \brief
398      * Returns the average histogram over all frames.
399      *
400      * Can be called already before the histogram is calculated to
401      * customize the way the average histogram is calculated.
402      *
403      * \see AbstractAverageHistogram
404      */
405     AbstractAverageHistogram& averager();
406
407     //! Returns bin properties for the histogram.
408     const AnalysisHistogramSettings& settings() const;
409
410     int frameCount() const override;
411
412     int flags() const override;
413
414     bool parallelDataStarted(AbstractAnalysisData* data, const AnalysisDataParallelOptions& options) override;
415     void frameStarted(const AnalysisDataFrameHeader& header) override;
416     void pointsAdded(const AnalysisDataPointSetRef& points) override;
417     void frameFinished(const AnalysisDataFrameHeader& header) override;
418     void frameFinishedSerial(int frameIndex) override;
419     void dataFinished() override;
420
421 private:
422     AnalysisDataFrameRef tryGetDataFrameInternal(int index) const override;
423     bool                 requestStorageInternal(int nframes) override;
424
425     class Impl;
426
427     PrivateImplPointer<Impl> impl_;
428
429     // Copy and assign disallowed by base.
430 };
431
432
433 /*! \brief
434  * Data module for per-frame weighted histograms.
435  *
436  * Output data contains the same number of frames and data sets as the input
437  * data.  Each frame contains the histogram(s) for the points in that frame,
438  * interpreted such that the first column passed to pointsAdded() determines
439  * the bin and the rest give weights to be added to that bin (input data should
440  * have at least two colums, and at least two columns should be added at the
441  * same time).
442  * Each input data set is processed independently into the corresponding output
443  * data set.
444  * All input columns for a data set are averaged into the same histogram.
445  * The number of columns for all data sets equals the number of bins in the
446  * histogram.
447  *
448  * The histograms are accumulated in double precision, even if the output data
449  * is in single precision.
450  *
451  * \inpublicapi
452  * \ingroup module_analysisdata
453  */
454 class AnalysisDataWeightedHistogramModule : public AbstractAnalysisData, public AnalysisDataModuleParallel
455 {
456 public:
457     //! \copydoc AnalysisDataSimpleHistogramModule::AnalysisDataSimpleHistogramModule()
458     AnalysisDataWeightedHistogramModule();
459     //! \copydoc AnalysisDataSimpleHistogramModule::AnalysisDataSimpleHistogramModule(const AnalysisHistogramSettings &)
460     explicit AnalysisDataWeightedHistogramModule(const AnalysisHistogramSettings& settings);
461     ~AnalysisDataWeightedHistogramModule() override;
462
463     //! \copydoc AnalysisDataSimpleHistogramModule::init()
464     void init(const AnalysisHistogramSettings& settings);
465
466     //! \copydoc AnalysisDataSimpleHistogramModule::averager()
467     AbstractAverageHistogram& averager();
468
469     //! \copydoc AnalysisDataSimpleHistogramModule::settings()
470     const AnalysisHistogramSettings& settings() const;
471
472     int frameCount() const override;
473
474     int flags() const override;
475
476     bool parallelDataStarted(AbstractAnalysisData* data, const AnalysisDataParallelOptions& options) override;
477     void frameStarted(const AnalysisDataFrameHeader& header) override;
478     void pointsAdded(const AnalysisDataPointSetRef& points) override;
479     void frameFinished(const AnalysisDataFrameHeader& header) override;
480     void frameFinishedSerial(int frameIndex) override;
481     void dataFinished() override;
482
483 private:
484     AnalysisDataFrameRef tryGetDataFrameInternal(int index) const override;
485     bool                 requestStorageInternal(int nframes) override;
486
487     class Impl;
488
489     PrivateImplPointer<Impl> impl_;
490
491     // Copy and assign disallowed by base.
492 };
493
494
495 /*! \brief
496  * Data module for bin averages.
497  *
498  * Output data contains one row for each bin; see AbstractAverageHistogram.
499  * Output data contains one column for each input data set.
500  * The value in a column is the average over all frames of that data set for
501  * that bin.
502  * The input data is interpreted such that the first column passed to
503  * pointsAdded() determines the bin and the rest give values to be added to
504  * that bin (input data should have at least two colums, and at least two
505  * columns should be added at the same time).
506  * All input columns for a data set are averaged into the same histogram.
507  *
508  * \inpublicapi
509  * \ingroup module_analysisdata
510  */
511 class AnalysisDataBinAverageModule : public AbstractAnalysisArrayData, public AnalysisDataModuleSerial
512 {
513 public:
514     //! \copydoc AnalysisDataSimpleHistogramModule::AnalysisDataSimpleHistogramModule()
515     AnalysisDataBinAverageModule();
516     //! \copydoc AnalysisDataSimpleHistogramModule::AnalysisDataSimpleHistogramModule(const AnalysisHistogramSettings &)
517     explicit AnalysisDataBinAverageModule(const AnalysisHistogramSettings& settings);
518     ~AnalysisDataBinAverageModule() override;
519
520     //! \copydoc AnalysisDataSimpleHistogramModule::init()
521     void init(const AnalysisHistogramSettings& settings);
522
523     //! \copydoc AnalysisDataSimpleHistogramModule::settings()
524     const AnalysisHistogramSettings& settings() const;
525
526     int flags() const override;
527
528     void dataStarted(AbstractAnalysisData* data) override;
529     void frameStarted(const AnalysisDataFrameHeader& header) override;
530     void pointsAdded(const AnalysisDataPointSetRef& points) override;
531     void frameFinished(const AnalysisDataFrameHeader& header) override;
532     void dataFinished() override;
533
534 private:
535     class Impl;
536
537     PrivateImplPointer<Impl> impl_;
538
539     // Copy and assign disallowed by base.
540 };
541
542 //! Smart pointer to manage an AnalysisDataSimpleHistogramModule object.
543 typedef std::shared_ptr<AnalysisDataSimpleHistogramModule> AnalysisDataSimpleHistogramModulePointer;
544 //! Smart pointer to manage an AnalysisDataWeightedHistogramModule object.
545 typedef std::shared_ptr<AnalysisDataWeightedHistogramModule> AnalysisDataWeightedHistogramModulePointer;
546 //! Smart pointer to manage an AnalysisDataBinAverageModule object.
547 typedef std::shared_ptr<AnalysisDataBinAverageModule> AnalysisDataBinAverageModulePointer;
548
549 } // namespace gmx
550
551 #endif