H2D.cpp
Go to the documentation of this file.
1 #ifdef __ICC
2 // disable icc remark #2259: non-pointer conversion from "X" to "Y" may lose significant bits
3 // TODO: To be removed, since it comes from ROOT TMathBase.h
4 #pragma warning(disable:2259)
5 #endif
6 #ifdef WIN32
7 // Disable warning
8 // warning C4996: 'sprintf': This function or variable may be unsafe.
9 // coming from TString.h
10 #pragma warning(disable:4996)
11 #endif
12 #include "H2D.h"
13 #include "H1D.h"
14 #include "P1D.h"
15 #include "TH1D.h"
16 #include "TProfile.h"
17 #include "GaudiKernel/ObjectFactory.h"
18 
19 std::pair<DataObject*,IHistogram2D*> Gaudi::createH2D(const std::string & title,int binsX,double iminX,double imaxX,int binsY,double iminY,double imaxY) {
20  auto p = new Histogram2D(new TH2D(title.c_str(),title.c_str(),binsX, iminX, imaxX, binsY, iminY, imaxY));
21  return { p, p };
22 }
23 
24 std::pair<DataObject*,IHistogram2D*> Gaudi::createH2D(const std::string & title, const Edges& eX, const Edges& eY) {
25  auto p = new Histogram2D(new TH2D(title.c_str(),title.c_str(),eX.size()-1,&eX.front(),eY.size()-1,&eY.front()));
26  return { p, p };
27 }
28 
29 std::pair<DataObject*,IHistogram2D*> Gaudi::createH2D(TH2D* rep) {
30  auto p = new Histogram2D(rep);
31  return { p, p };
32 }
33 
34 std::pair<DataObject*,IHistogram2D*> Gaudi::createH2D(const IHistogram2D& hist) {
35  TH2D *h = getRepresentation<AIDA::IHistogram2D,TH2D>(hist);
36  Histogram2D *n = h ? new Histogram2D(new TH2D(*h)) : nullptr;
37  return { n, n };
38 }
39 
40 std::pair<DataObject*,IHistogram1D*>
41 Gaudi::slice1DX(const std::string& nam,const IHistogram2D& hist,int first, int last) {
42  TH2 *r = getRepresentation<IHistogram2D,TH2>(hist);
43  TH1D *t = r ? r->ProjectionX("_px",first,last,"e") : nullptr;
44  if ( t ) t->SetName(nam.c_str());
45  Histogram1D* p = ( t ? new Histogram1D(t) : nullptr );
46  return { p, p };
47 }
48 
49 std::pair<DataObject*,IHistogram1D*>
50 Gaudi::slice1DY(const std::string& nam,const IHistogram2D& hist,int first, int last) {
51  TH2 *r = getRepresentation<IHistogram2D,TH2>(hist);
52  TH1D *t = r ? r->ProjectionY("_py",first,last,"e") : nullptr;
53  if ( t ) t->SetName(nam.c_str());
54  Histogram1D* p = ( t ? new Histogram1D(t) : nullptr );
55  return { p, p };
56 }
57 
58 std::pair<DataObject*,IHistogram1D*>
59 Gaudi::project1DY(const std::string& nam,const IHistogram2D& hist,int first,int last) {
60  TH2 *r = getRepresentation<IHistogram2D,TH2>(hist);
61  TH1D *t = r ? r->ProjectionY("_px",first,last,"e") : nullptr;
62  if ( t ) t->SetName(nam.c_str());
63  Histogram1D* p = ( t ? new Histogram1D(t) : nullptr );
64  return { p, p };
65 }
66 
67 std::pair<DataObject*,IProfile1D*>
68 Gaudi::profile1DX(const std::string& nam,const IHistogram2D& hist,int first,int last) {
69  TH2 *r = Gaudi::getRepresentation<IHistogram2D,TH2>(hist);
70  TProfile *t = r ? r->ProfileX("_pfx",first,last,"e") : nullptr;
71  if ( t ) t->SetName(nam.c_str());
72  Profile1D* p = ( t ? new Profile1D(t) : nullptr );
73  return { p, p };
74 }
75 
76 std::pair<DataObject*,IProfile1D*>
77 Gaudi::profile1DY(const std::string& nam, const IHistogram2D& hist,int first, int last) {
78  TH2 *r = getRepresentation<IHistogram2D,TH2>(hist);
79  TProfile *t = r ? r->ProfileY("_pfx",first,last,"e") : nullptr;
80  if ( t ) t->SetName(nam.c_str());
81  Profile1D* p = ( t ? new Profile1D(t) : nullptr );
82  return { p, p };
83 }
84 
85 
86 namespace Gaudi {
87  template <>
88  void * Generic2D<IHistogram2D,TH2D>::cast(const std::string & className) const {
89  if (className == "AIDA::IHistogram2D")
90  return (IHistogram2D*)this;
91  else if (className == "AIDA::IHistogram")
92  return (IHistogram*)this;
93  return nullptr;
94  }
95 
96  template <>
97  int Generic2D<IHistogram2D,TH2D>::binEntries(int indexX,int indexY) const {
98  if (binHeight(indexX, indexY)<=0) return 0;
99  double xx = binHeight(indexX, indexY)/binError(indexX, indexY);
100  return int(xx*xx+0.5);
101  }
102 
103  template <>
105  TH2D* imp = dynamic_cast<TH2D*>(rep);
106  if ( !imp ) throw std::runtime_error("Cannot adopt native histogram representation.");
107  m_rep.reset( imp );
108  m_xAxis.initialize(m_rep->GetXaxis(),true);
109  m_yAxis.initialize(m_rep->GetYaxis(),true);
110  const TArrayD* a = m_rep->GetSumw2();
111  if ( !a || (a && a->GetSize()==0) ) m_rep->Sumw2();
112  setTitle(m_rep->GetTitle());
113  }
114 }
115 
117  : Base( new TH2D() )
118 {
119  m_rep->Sumw2();
120  m_sumwx = m_sumwy = 0;
121  setTitle("");
122  m_rep->SetDirectory(nullptr);
123 }
124 
126 {
127  adoptRepresentation(rep);
128  m_sumwx = m_sumwy = 0;
129  m_rep->SetDirectory(nullptr);
130 }
131 
132 bool Gaudi::Histogram2D::setBinContents( int i,int j,int entries,double height,double error,double centreX,double centreY) {
133  m_rep->SetBinContent(rIndexX(i), rIndexY(j), height);
134  m_rep->SetBinError(rIndexX(i), rIndexY(j), error);
135  // accumulate sumwx for in range bins
136  if (i >=0 && j >= 0) {
137  m_sumwx += centreX*height;
138  m_sumwy += centreY*height;
139  }
140  m_sumEntries += entries;
141  return true;
142 }
143 
145  m_sumwx = 0;
146  m_sumwy = 0;
147  return Base::reset();
148 }
149 
150 #ifdef __ICC
151 // disable icc remark #1572: floating-point equality and inequality comparisons are unreliable
152 // The comparison are meant
153 #pragma warning(push)
154 #pragma warning(disable:1572)
155 #endif
156 bool Gaudi::Histogram2D::fill ( double x,double y,double weight) {
157  (weight==1.) ? m_rep->Fill(x,y) : m_rep->Fill(x,y,weight );
158  return true;
159 }
160 
161 bool Gaudi::Histogram2D::setRms(double rmsX,double rmsY) {
162  m_rep->SetEntries(m_sumEntries);
163  std::vector<double> stat(11);
164  stat[0] = sumBinHeights();
165  stat[1] = 0;
166  if(equivalentBinEntries() != 0)
167  stat[1] = (sumBinHeights() * sumBinHeights()) / equivalentBinEntries();
168  stat[2] = m_sumwx;
169  double meanX = 0;
170  if(sumBinHeights() != 0) meanX = m_sumwx/ sumBinHeights();
171  stat[3] = (meanX*meanX + rmsX*rmsX) * sumBinHeights();
172  stat[4] = m_sumwy;
173  double meanY = 0;
174  if(sumBinHeights() != 0) meanY = m_sumwy/ sumBinHeights();
175  stat[5] = (meanY*meanY + rmsY*rmsY) * sumBinHeights();
176  stat[6] = 0;
177  m_rep->PutStats(&stat.front());
178  return true;
179 }
180 
181 void Gaudi::Histogram2D::copyFromAida(const IHistogram2D& h) {
182  // implement here the copy
183  const char* tit = h.title().c_str();
184  if (h.xAxis().isFixedBinning() && h.yAxis().isFixedBinning() )
185  m_rep.reset( new TH2D(tit,tit,
186  h.xAxis().bins(),h.xAxis().lowerEdge(),h.xAxis().upperEdge(),
187  h.yAxis().bins(),h.yAxis().lowerEdge(),h.yAxis().upperEdge() ) );
188  else {
189  Edges eX, eY;
190  for (int i =0; i < h.xAxis().bins(); ++i)
191  eX.push_back(h.xAxis().binLowerEdge(i));
192  // add also upperedges at the end
193  eX.push_back(h.xAxis().upperEdge() );
194  for (int i =0; i < h.yAxis().bins(); ++i)
195  eY.push_back(h.yAxis().binLowerEdge(i));
196  // add also upperedges at the end
197  eY.push_back(h.yAxis().upperEdge() );
198  m_rep.reset( new TH2D(tit,tit,eX.size()-1,&eX.front(),eY.size()-1,&eY.front()) );
199  }
200  m_xAxis.initialize(m_rep->GetXaxis(),true);
201  m_yAxis.initialize(m_rep->GetYaxis(),true);
202  m_rep->Sumw2();
203  m_sumEntries = 0;
204  m_sumwx = 0;
205  m_sumwy = 0;
206  // statistics
207  double sumw = h.sumBinHeights();
208  double sumw2 = 0;
209  if (h.equivalentBinEntries() != 0)
210  sumw2 = ( sumw * sumw ) /h.equivalentBinEntries();
211  double sumwx = h.meanX()*h.sumBinHeights();
212  double sumwx2 = (h.meanX()*h.meanX() + h.rmsX()*h.rmsX() )*h.sumBinHeights();
213  double sumwy = h.meanY()*h.sumBinHeights();
214  double sumwy2 = (h.meanY()*h.meanY() + h.rmsY()*h.rmsY() )*h.sumBinHeights();
215  double sumwxy = 0;
216 
217  // copy the contents in (AIDA underflow/overflow are -2,-1)
218  for (int i=-2; i < xAxis().bins(); ++i) {
219  for (int j=-2; j < yAxis().bins(); ++j) {
220  // root binning starts from one !
221  m_rep->SetBinContent(rIndexX(i), rIndexY(j), h.binHeight(i,j) );
222  m_rep->SetBinError(rIndexX(i), rIndexY(j), h.binError(i,j) );
223  // calculate statistics
224  if ( i >= 0 && j >= 0) {
225  sumwxy += h.binHeight(i,j)*h.binMeanX(i,j)*h.binMeanY(i,j);
226  }
227  }
228  }
229  // need to do set entries after setting contents otherwise root will recalculate them
230  // taking into account how many time SetBinContents() has been called
231  m_rep->SetEntries(h.allEntries());
232  // fill stat vector
233  std::array<double,11> stat = {{ sumw, sumw2,
234  sumwx, sumwx2,
235  sumwy, sumwy2,
236  sumwxy }};
237  m_rep->PutStats(stat.data());
238 }
239 
242 
243 #ifdef __ICC
244 // re-enable icc remark #1572
245 #pragma warning(pop)
246 #endif
std::pair< DataObject *, AIDA::IHistogram1D * > slice1DY(const std::string &name, const AIDA::IHistogram2D &h, int firstbin, int lastbin)
Create 1D slice from 2D histogram.
bool setRms(double rmsX, double rmsY)
Sets the rms of the histogram.
Definition: H2D.cpp:161
std::pair< DataObject *, AIDA::IHistogram1D * > project1DY(const std::string &name, const AIDA::IHistogram2D &h, int firstbin, int lastbin)
Create 1D projection in Y from 2D histogram.
virtual bool setBinContents(int binIndexX, int binIndexY, int entries, double height, double error, double centreX, double centreY)
Fast filling method for a given bin. It can be also the over/underflow bin.
Definition: H2D.cpp:132
std::vector< double > Edges
Definition: GaudiPI.h:19
int binEntries(int indexX, int indexY) const override
The number of entries (ie the number of times fill was called for this bin).
void adoptRepresentation(TObject *rep) override
Adopt ROOT histogram representation.
Histogram2D()
Standard Constructor.
Definition: H2D.cpp:116
void * cast(const std::string &className) const
Introspection method.
AIDA implementation for 2 D histograms using ROOT THD2.
Definition: H2D.h:20
bool reset()
Definition: H2D.cpp:144
std::pair< DataObject *, AIDA::IHistogram2D * > createH2D(const AIDA::IHistogram2D &hist)
Copy constructor.
std::pair< DataObject *, AIDA::IHistogram1D * > slice1DX(const std::string &name, const AIDA::IHistogram2D &h, int firstbin, int lastbin)
Create 1D slice from 2D histogram.
std::pair< DataObject *, AIDA::IProfile1D * > profile1DX(const std::string &name, const AIDA::IHistogram2D &h, int firstbin, int lastbin)
Create 1D profile in X from 2D histogram.
std::pair< DataObject *, AIDA::IProfile1D * > profile1DY(const std::string &name, const AIDA::IHistogram2D &h, int firstbin, int lastbin)
Create 1D profile in Y from 2D histogram.
bool setTitle(const std::string &title) override
Set the title of the object.
Definition: Generic2D.h:158
Gaudi::Histogram2D H2D
Definition: H2D.cpp:240
bool fill(double x, double y, double weight=1.)
Fill the Histogram2D with a value and the.
Definition: H2D.cpp:156
void copyFromAida(const IHistogram2D &h)
Create new histogram from any AIDA based histogram.
Definition: H2D.cpp:181
std::unique_ptr< IMPLEMENTATION > m_rep
Reference to underlying implementation.
Definition: Generic2D.h:150
double m_sumwy
Definition: H2D.h:46
double m_sumwx
Definition: H2D.h:45
#define DECLARE_DATAOBJECT_FACTORY(x)
Definition: ObjectFactory.h:21
list i
Definition: ana.py:128
Helper functions to set/get the application return code.
Definition: __init__.py:1
Common AIDA implementation stuff for histograms and profiles using ROOT implementations.
Definition: Generic2D.h:35