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