Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ConvertToTH1.cxx
Go to the documentation of this file.
1/// \file
2/// \warning This is part of the %ROOT 7 prototype! It will change without notice. It might trigger earthquakes.
3/// Feedback is welcome!
4
7#include <ROOT/RBinIndex.hxx>
8#include <ROOT/RHist.hxx>
10
11#include <TH1.h>
12
13#include <cassert>
14#include <cstddef>
15#include <memory>
16#include <stdexcept>
17#include <type_traits>
18#include <variant>
19#include <vector>
20
21using namespace ROOT::Experimental;
22
23namespace {
24template <typename Hist, typename T>
25std::unique_ptr<Hist> ConvertToTH1Impl(const RHistEngine<T> &engine)
26{
27 if (engine.GetNDimensions() != 1) {
28 throw std::invalid_argument("TH1 requires one dimension");
29 }
30
31 auto ret = std::make_unique<Hist>();
32 ret->SetDirectory(nullptr);
33
34 ROOT::Experimental::Hist::Internal::ConvertAxis(*ret->GetXaxis(), engine.GetAxes()[0]);
35 ret->SetBinsLength();
36
37 Double_t *sumw2 = nullptr;
38 auto copyBinContent = [&ret, &engine, &sumw2](Int_t i, RBinIndex index) {
39 if constexpr (std::is_same_v<T, RBinWithError>) {
40 if (sumw2 == nullptr) {
41 ret->Sumw2();
42 sumw2 = ret->GetSumw2()->GetArray();
43 }
44 const RBinWithError &c = engine.GetBinContent(index);
45 ret->GetArray()[i] = c.fSum;
46 sumw2[i] = c.fSum2;
47 } else {
48 (void)sumw2;
49 ret->GetArray()[i] = engine.GetBinContent(index);
50 }
51 };
52
53 // Copy the bin contents, accounting for TH1 numbering conventions.
54 const auto &axis = engine.GetAxes()[0];
55 for (auto index : axis.GetFullRange()) {
56 if (index.IsUnderflow()) {
58 } else if (index.IsOverflow()) {
59 copyBinContent(axis.GetNNormalBins() + 1, index);
60 } else {
61 assert(index.IsNormal());
62 copyBinContent(index.GetIndex() + 1, index);
63 }
64 }
65
66 return ret;
67}
68
69template <typename Hist>
70void ConvertGlobalStatistics(Hist &h, const RHistStats &stats)
71{
72 if (stats.IsTainted()) {
73 return;
74 }
75
76 h.SetEntries(stats.GetNEntries());
77
78 Double_t hStats[4] = {
79 stats.GetSumW(),
80 stats.GetSumW2(),
81 0,
82 0,
83 };
84 if (stats.IsEnabled(0)) {
85 hStats[2] = stats.GetDimensionStats(0).fSumWX;
86 hStats[3] = stats.GetDimensionStats(0).fSumWX2;
87 }
88 h.PutStats(hStats);
89}
90} // namespace
91
92namespace ROOT {
93namespace Experimental {
94namespace Hist {
95
96std::unique_ptr<TH1C> ConvertToTH1C(const RHistEngine<char> &engine)
97{
98 return ConvertToTH1Impl<TH1C>(engine);
99}
100
101std::unique_ptr<TH1S> ConvertToTH1S(const RHistEngine<short> &engine)
102{
103 return ConvertToTH1Impl<TH1S>(engine);
104}
105
106std::unique_ptr<TH1I> ConvertToTH1I(const RHistEngine<int> &engine)
107{
108 return ConvertToTH1Impl<TH1I>(engine);
109}
110
111std::unique_ptr<TH1L> ConvertToTH1L(const RHistEngine<long> &engine)
112{
113 return ConvertToTH1Impl<TH1L>(engine);
114}
115
116std::unique_ptr<TH1L> ConvertToTH1L(const RHistEngine<long long> &engine)
117{
118 return ConvertToTH1Impl<TH1L>(engine);
119}
120
121std::unique_ptr<TH1F> ConvertToTH1F(const RHistEngine<float> &engine)
122{
123 return ConvertToTH1Impl<TH1F>(engine);
124}
125
126std::unique_ptr<TH1D> ConvertToTH1D(const RHistEngine<double> &engine)
127{
128 return ConvertToTH1Impl<TH1D>(engine);
129}
130
131std::unique_ptr<TH1D> ConvertToTH1D(const RHistEngine<RBinWithError> &engine)
132{
133 return ConvertToTH1Impl<TH1D>(engine);
134}
135
136std::unique_ptr<TH1C> ConvertToTH1C(const RHist<char> &hist)
137{
138 auto ret = ConvertToTH1C(hist.GetEngine());
139 ConvertGlobalStatistics(*ret, hist.GetStats());
140 return ret;
141}
142
143std::unique_ptr<TH1S> ConvertToTH1S(const RHist<short> &hist)
144{
145 auto ret = ConvertToTH1S(hist.GetEngine());
146 ConvertGlobalStatistics(*ret, hist.GetStats());
147 return ret;
148}
149
150std::unique_ptr<TH1I> ConvertToTH1I(const RHist<int> &hist)
151{
152 auto ret = ConvertToTH1I(hist.GetEngine());
153 ConvertGlobalStatistics(*ret, hist.GetStats());
154 return ret;
155}
156
157std::unique_ptr<TH1L> ConvertToTH1L(const RHist<long> &hist)
158{
159 auto ret = ConvertToTH1L(hist.GetEngine());
160 ConvertGlobalStatistics(*ret, hist.GetStats());
161 return ret;
162}
163
164std::unique_ptr<TH1L> ConvertToTH1L(const RHist<long long> &hist)
165{
166 auto ret = ConvertToTH1L(hist.GetEngine());
167 ConvertGlobalStatistics(*ret, hist.GetStats());
168 return ret;
169}
170
171std::unique_ptr<TH1F> ConvertToTH1F(const RHist<float> &hist)
172{
173 auto ret = ConvertToTH1F(hist.GetEngine());
174 ConvertGlobalStatistics(*ret, hist.GetStats());
175 return ret;
176}
177
178std::unique_ptr<TH1D> ConvertToTH1D(const RHist<double> &hist)
179{
180 auto ret = ConvertToTH1D(hist.GetEngine());
181 ConvertGlobalStatistics(*ret, hist.GetStats());
182 return ret;
183}
184
185std::unique_ptr<TH1D> ConvertToTH1D(const RHist<RBinWithError> &hist)
186{
187 auto ret = ConvertToTH1D(hist.GetEngine());
188 ConvertGlobalStatistics(*ret, hist.GetStats());
189 return ret;
190}
191
192} // namespace Hist
193} // namespace Experimental
194} // namespace ROOT
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
A bin index with special values for underflow and overflow bins.
Definition RBinIndex.hxx:23
Histogram statistics of unbinned values.
const RDimensionStats & GetDimensionStats(std::size_t dim=0) const
Get the statistics object for one dimension.
std::uint64_t GetNEntries() const
bool IsEnabled(std::size_t dim) const
void ConvertAxis(TAxis &dst, const RAxisVariant &src)
Convert a single axis object to TAxis.
std::unique_ptr< TH1I > ConvertToTH1I(const RHistEngine< int > &engine)
Convert a one-dimensional histogram to TH1I.
std::unique_ptr< TH1S > ConvertToTH1S(const RHistEngine< short > &engine)
Convert a one-dimensional histogram to TH1S.
std::unique_ptr< TH1L > ConvertToTH1L(const RHistEngine< long > &engine)
Convert a one-dimensional histogram to TH1L.
std::unique_ptr< TH1F > ConvertToTH1F(const RHistEngine< float > &engine)
Convert a one-dimensional histogram to TH1F.
std::unique_ptr< TH1C > ConvertToTH1C(const RHistEngine< char > &engine)
Convert a one-dimensional histogram to TH1C.
std::unique_ptr< TH1D > ConvertToTH1D(const RHistEngine< double > &engine)
Convert a one-dimensional histogram to TH1D.
Namespace for ROOT features in testing.
Definition TROOT.h:100