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