Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
THnChain.cxx
Go to the documentation of this file.
1// @(#)root/hist:$Id$
2// Author: Benjamin Bannier, August 2016
3
4/*************************************************************************
5 * Copyright (C) 1995-2016, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#include "THnChain.h"
13
14#include "TArray.h"
15#include "TAxis.h"
16#include "TDirectory.h"
17#include "TFile.h"
18#include "TH1.h"
19#include "TH2.h"
20#include "TH3.h"
21#include "THnBase.h"
22#include "TMath.h"
23
24/// Add a new file to this chain.
25///
26/// \param fileName path of the file to add
27void THnChain::AddFile(const char* fileName)
28{
29 fFiles.emplace_back(fileName);
30
31 // Initialize axes from first seen instance.
32 if (fAxes.empty()) {
33 THnBase* hs = ReadHistogram(fileName);
34
35 if (hs) {
36 const Int_t naxes = hs->GetNdimensions();
37 for (Int_t i = 0; i < naxes; ++i) {
38 fAxes.push_back(hs->GetAxis(i));
39 }
40 } else {
41 Warning("AddFile",
42 "Could not find histogram %s in file %s",
43 fName.c_str(),
44 fileName);
45 }
46 }
47}
48
49/// Get an axis from the histogram.
50///
51/// \param i index of the axis to retrieve
52///
53/// This function requires that a file containing the histogram was
54/// already added with `AddFile`.
55///
56/// Properties set on the axis returned by `GetAxis` are propagated to all
57/// histograms in the chain, so it can e.g. be used to set ranges for
58/// projections.
60{
61 if (i < 0 || i >= static_cast<Int_t>(fAxes.size())) {
62 return nullptr;
63 }
64
65 return fAxes[i];
66}
67
68/// Projects all histograms in the chain.
69///
70/// See `THnBase::Projection` for parameters and their semantics.
72{
73 if (ndim <= 0) {
74 return nullptr;
75 }
76
77 TObject* h_merged = nullptr;
78 for (const auto& file : fFiles) {
79 THnBase* hs = ReadHistogram(file.c_str());
80
81 if (!hs) {
82 Warning("ProjectionAny",
83 "Could not find histogram %s in file %s",
84 fName.c_str(),
85 file.c_str());
86
87 continue;
88 }
89
90 if (!CheckConsistency(*hs, fAxes)) {
91 Warning("ProjectionAny",
92 "Histogram %s from file %s is inconsistent with the histogram from file %s",
93 fName.c_str(),
94 file.c_str(),
95 fFiles[0].c_str());
96
97 continue;
98 }
99
100 SetupAxes(*hs);
101
102 // Perform projection.
103 TObject* h = nullptr;
104
105 if (ndim == 1) {
106 h = hs->Projection(dim[0], option);
107 } else if (ndim == 2) {
108 h = hs->Projection(dim[0], dim[1], option);
109 } else if (ndim == 3) {
110 h = hs->Projection(dim[0], dim[1], dim[2], option);
111 } else {
112 h = hs->ProjectionND(ndim, dim, option);
113 }
114
115 delete hs;
116
117 // Add this histogram.
118 if (h_merged) {
119 if (ndim < 3) {
120 static_cast<TH1*>(h_merged)->Add(static_cast<TH1*>(h));
121 } else {
122 static_cast<THnBase*>(h_merged)->Add(static_cast<THnBase*>(h));
123 }
124
125 delete h;
126 } else {
127 h_merged = h;
128 }
129 }
130
131 return h_merged;
132}
133
134/// Retrieve a histogram from a file.
135///
136/// \param fileName path of the file to read.
137THnBase* THnChain::ReadHistogram(const char* fileName) const
138{
140
141 TFile* f = TFile::Open(fileName);
142
143 if (!f) {
144 return nullptr;
145 }
146
147 THnBase* hs = nullptr;
148 f->GetObject(fName.c_str(), hs);
149 delete f;
150
151 return hs;
152}
153
154/// Copy the properties of all axes to a histogram.
155///
156/// \param hs histogram whose axes should be updated
158{
159 const Int_t naxes = fAxes.size();
160 for (Int_t i = 0; i < naxes; ++i) {
161 const TAxis* ax_ref = fAxes[i];
162 TAxis* ax = hs.GetAxis(i);
163 ax_ref->Copy(*ax);
164 }
165}
166
167/// Ensure a histogram has axes similar to the ones we expect.
168///
169/// \param h histogram to verify
170/// \param axes expected set of axes
171bool THnChain::CheckConsistency(const THnBase& h, const std::vector<TAxis*>& axes)
172{
173 // We would prefer to directly use `TH1::CheckEqualAxes` here;
174 // however it is protected so we inherit the parts we care about.
175 // FIXME(bbannier): It appears that functionality like `TH1::CheckEqualAxes` could
176 // just as well live in `TAxis` so that anyone using axes could make use of it.
177 const Int_t naxes = h.GetNdimensions();
178 const Int_t naxes2 = axes.size();
179
180 if (naxes != naxes2) {
181 return false;
182 }
183
184 for (Int_t i = 0; i < naxes; ++i) {
185 const TAxis* ax1 = h.GetAxis(i);
186 const TAxis* ax2 = axes[i];
187
188 if (ax1->GetNbins() != ax2->GetNbins()) {
189 return false;
190 }
191
192 // Copied from `TH1::CheckAxisLimits.`
193 if (!TMath::AreEqualRel(ax1->GetXmin(), ax2->GetXmin(), 1.E-12) ||
194 !TMath::AreEqualRel(ax1->GetXmax(), ax2->GetXmax(), 1.E-12)) {
195 return false;
196 }
197
198 // Copied from `TH1::CheckBinLimits`.
199 const TArrayD* h1Array = ax1->GetXbins();
200 const TArrayD* h2Array = ax2->GetXbins();
201 Int_t fN = h1Array->fN;
202 if (fN != 0) {
203 if (h2Array->fN != fN) {
204 return false;
205 } else {
206 for (int ibin = 0; ibin < fN; ++ibin) {
207 if (!TMath::AreEqualRel(h1Array->GetAt(ibin), h2Array->GetAt(ibin), 1E-10)) {
208 return false;
209 }
210 }
211 }
212 }
213
214 // We ignore checking for equal bin labels here. A check
215 // for that is implemented in `TH1::CheckBinLabels`.
216 }
217
218 return true;
219}
220
221/// See `THnBase::Projection` for the intended behavior.
223{
224 // Forwards to `THnBase::Projection()`.
225 Int_t dim[1] = {xDim};
226 return static_cast<TH1*>(ProjectionAny(1, dim, option));
227}
228
229/// See `THnBase::Projection` for the intended behavior.
231{
232 // Forwards to `THnBase::Projection()`.
233 Int_t dim[2] = {xDim, yDim};
234 return static_cast<TH2*>(ProjectionAny(2, dim, option));
235}
236
237/// See `THnBase::Projection` for the intended behavior.
239{
240 // Forwards to `THnBase::Projection()`.
241 Int_t dim[3] = {xDim, yDim, zDim};
242 return static_cast<TH3*>(ProjectionAny(3, dim, option));
243}
244
245/// See `THnBase::Projection` for the intended behavior.
247{
248 // Forwards to `THnBase::ProjectionND()`.
249 return static_cast<THnBase*>(ProjectionAny(ndim, dim, option));
250}
#define f(i)
Definition RSha256.hxx:104
#define h(i)
Definition RSha256.hxx:106
const char Option_t
Definition RtypesCore.h:66
#define gDirectory
Definition TDirectory.h:384
Option_t Option_t option
Array of doubles (64 bits per element).
Definition TArrayD.h:27
Double_t GetAt(Int_t i) const override
Definition TArrayD.h:45
Int_t fN
Definition TArray.h:38
Class to manage histogram axis.
Definition TAxis.h:31
const TArrayD * GetXbins() const
Definition TAxis.h:136
void Copy(TObject &axis) const override
Copy axis structure to another axis.
Definition TAxis.cxx:216
Double_t GetXmax() const
Definition TAxis.h:140
Double_t GetXmin() const
Definition TAxis.h:139
Int_t GetNbins() const
Definition TAxis.h:125
TDirectory::TContext keeps track and restore the current directory.
Definition TDirectory.h:89
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:53
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:4089
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
Service class for 2-D histogram classes.
Definition TH2.h:30
The 3-D histogram classes derived from the 1-D histogram classes.
Definition TH3.h:31
Multidimensional histogram base.
Definition THnBase.h:43
TH1D * Projection(Int_t xDim, Option_t *option="") const
Project all bins into a 1-dimensional histogram, keeping only axis "xDim".
Definition THnBase.h:226
Int_t GetNdimensions() const
Definition THnBase.h:140
THnBase * ProjectionND(Int_t ndim, const Int_t *dim, Option_t *option="") const
Definition THnBase.h:257
TAxis * GetAxis(Int_t dim) const
Definition THnBase.h:130
TAxis * GetAxis(Int_t i) const
Get an axis from the histogram.
Definition THnChain.cxx:59
std::vector< TAxis * > fAxes
the list of histogram axes
Definition THnChain.h:73
THnBase * ReadHistogram(const char *fileName) const
Retrieve a histogram from a file.
Definition THnChain.cxx:137
void SetupAxes(THnBase &hs) const
Copy the properties of all axes to a histogram.
Definition THnChain.cxx:157
THnBase * ProjectionND(Int_t ndim, const Int_t *dim, Option_t *option="") const
See THnBase::Projection for the intended behavior.
Definition THnChain.cxx:246
TObject * ProjectionAny(Int_t ndim, const Int_t *dim, Option_t *option="") const
Projects all histograms in the chain.
Definition THnChain.cxx:71
static bool CheckConsistency(const THnBase &h, const std::vector< TAxis * > &axes)
Ensure a histogram has axes similar to the ones we expect.
Definition THnChain.cxx:171
void AddFile(const char *fileName)
Add a new file to this chain.
Definition THnChain.cxx:27
TH1 * Projection(Int_t xDim, Option_t *option="") const
See THnBase::Projection for the intended behavior.
Definition THnChain.cxx:222
std::string fName
name of the histogram
Definition THnChain.h:70
std::vector< std::string > fFiles
a list of files to extract the histogram from
Definition THnChain.h:72
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:991
Bool_t AreEqualRel(Double_t af, Double_t bf, Double_t relPrec)
Comparing floating points.
Definition TMath.h:426
TMatrixT< Element > & Add(TMatrixT< Element > &target, Element scalar, const TMatrixT< Element > &source)
Modify addition: target += scalar * source.