Gaudi Framework, version v23r6

Home   Generated: Wed Jan 30 2013
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
H2D.cpp
Go to the documentation of this file.
1 // $Id: H2D.cpp,v 1.16 2007/07/16 13:36:17 hmd Exp $
2 #ifdef __ICC
3 // disable icc remark #2259: non-pointer conversion from "X" to "Y" may lose significant bits
4 // TODO: To be removed, since it comes from ROOT TMathBase.h
5 #pragma warning(disable:2259)
6 #endif
7 #ifdef WIN32
8 // Disable warning
9 // warning C4996: 'sprintf': This function or variable may be unsafe.
10 // coming from TString.h
11 #pragma warning(disable:4996)
12 #endif
13 #include "H2D.h"
14 #include "H1D.h"
15 #include "P1D.h"
16 #include "TH1D.h"
17 #include "TProfile.h"
19 
20 std::pair<DataObject*,IHistogram2D*> Gaudi::createH2D(const std::string & title,int binsX,double iminX,double imaxX,int binsY,double iminY,double imaxY) {
21  Histogram2D* p = new Histogram2D(new TH2D(title.c_str(),title.c_str(),binsX, iminX, imaxX, binsY, iminY, imaxY));
23 }
24 
25 std::pair<DataObject*,IHistogram2D*> Gaudi::createH2D(const std::string & title, const Edges& eX, const Edges& eY) {
26  Histogram2D* p = new Histogram2D(new TH2D(title.c_str(),title.c_str(),eX.size()-1,&eX.front(),eY.size()-1,&eY.front()));
28 }
29 
31  Histogram2D* p = new Histogram2D(rep);
33 }
34 
36  TH2D *h = getRepresentation<AIDA::IHistogram2D,TH2D>(hist);
37  Histogram2D *n = h ? new Histogram2D(new TH2D(*h)) : 0;
39 }
40 
42 Gaudi::slice1DX(const std::string& nam,const IHistogram2D& hist,int first, int last) {
43  TH2 *r = getRepresentation<IHistogram2D,TH2>(hist);
44  TH1D *t = r ? r->ProjectionX("_px",first,last,"e") : 0;
45  if ( t ) t->SetName(nam.c_str());
46  Histogram1D* p = t ? new Histogram1D(t) : 0;
48 }
49 
51 Gaudi::slice1DY(const std::string& nam,const IHistogram2D& hist,int first, int last) {
52  TH2 *r = getRepresentation<IHistogram2D,TH2>(hist);
53  TH1D *t = r ? r->ProjectionY("_py",first,last,"e") : 0;
54  if ( t ) t->SetName(nam.c_str());
55  Histogram1D* p = t ? new Histogram1D(t) : 0;
57 }
58 
60 Gaudi::project1DY(const std::string& nam,const IHistogram2D& hist,int first,int last) {
61  TH2 *r = getRepresentation<IHistogram2D,TH2>(hist);
62  TH1D *t = r ? r->ProjectionY("_px",first,last,"e") : 0;
63  if ( t ) t->SetName(nam.c_str());
64  Histogram1D* p = t ? new Histogram1D(t) : 0;
66 }
67 
69 Gaudi::profile1DX(const std::string& nam,const IHistogram2D& hist,int first,int last) {
70  TH2 *r = Gaudi::getRepresentation<IHistogram2D,TH2>(hist);
71  TProfile *t = r ? r->ProfileX("_pfx",first,last,"e") : 0;
72  if ( t ) t->SetName(nam.c_str());
73  Profile1D* p = t ? new Profile1D(t) : 0;
75 }
76 
78 Gaudi::profile1DY(const std::string& nam, const IHistogram2D& hist,int first, int last) {
79  TH2 *r = getRepresentation<IHistogram2D,TH2>(hist);
80  TProfile *t = r ? r->ProfileY("_pfx",first,last,"e") : 0;
81  if ( t ) t->SetName(nam.c_str());
82  Profile1D* p = t ? new Profile1D(t) : 0;
84 }
85 
86 
87 namespace Gaudi {
88  template <>
89  void * Generic2D<IHistogram2D,TH2D>::cast(const std::string & className) const {
90  if (className == "AIDA::IHistogram2D")
91  return (IHistogram2D*)this;
92  else if (className == "AIDA::IHistogram")
93  return (IHistogram*)this;
94  return 0;
95  }
96 
97  template <>
98  int Generic2D<IHistogram2D,TH2D>::binEntries(int indexX,int indexY) const {
99  if (binHeight(indexX, indexY)<=0) return 0;
100  double xx = binHeight(indexX, indexY)/binError(indexX, indexY);
101  return int(xx*xx+0.5);
102  }
103 
104  template <>
106  TH2D* imp = dynamic_cast<TH2D*>(rep);
107  if ( imp ) {
108  if ( m_rep ) delete m_rep;
109  m_rep = imp;
110  m_xAxis.initialize(m_rep->GetXaxis(),true);
111  m_yAxis.initialize(m_rep->GetYaxis(),true);
112  const TArrayD* a = m_rep->GetSumw2();
113  if ( 0 == a || (a && a->GetSize()==0) ) m_rep->Sumw2();
114  setTitle(m_rep->GetTitle());
115  return;
116  }
117  throw std::runtime_error("Cannot adopt native histogram representation.");
118  }
119 }
120 
122  m_rep = new TH2D();
123  m_rep->Sumw2();
124  m_sumEntries = 0;
125  m_sumwx = m_sumwy = 0;
126  setTitle("");
127  m_rep->SetDirectory(0);
128 }
129 
131  m_rep = 0;
132  adoptRepresentation(rep);
133  m_sumEntries = 0;
134  m_sumwx = m_sumwy = 0;
135  m_rep->SetDirectory(0);
136 }
137 
138 bool Gaudi::Histogram2D::setBinContents( int i,int j,int entries,double height,double error,double centreX,double centreY) {
139  m_rep->SetBinContent(rIndexX(i), rIndexY(j), height);
140  m_rep->SetBinError(rIndexX(i), rIndexY(j), error);
141  // accumulate sumwx for in range bins
142  if (i >=0 && j >= 0) {
143  m_sumwx += centreX*height;
144  m_sumwy += centreY*height;
145  }
146  m_sumEntries += entries;
147  return true;
148 }
149 
151  m_sumwx = 0;
152  m_sumwy = 0;
153  return Base::reset();
154 }
155 
156 #ifdef __ICC
157 // disable icc remark #1572: floating-point equality and inequality comparisons are unreliable
158 // The comparison are meant
159 #pragma warning(push)
160 #pragma warning(disable:1572)
161 #endif
162 bool Gaudi::Histogram2D::fill ( double x,double y,double weight) {
163  (weight==1.) ? m_rep->Fill(x,y) : m_rep->Fill(x,y,weight );
164  return true;
165 }
166 
167 bool Gaudi::Histogram2D::setRms(double rmsX,double rmsY) {
168  m_rep->SetEntries(m_sumEntries);
169  std::vector<double> stat(11);
170  stat[0] = sumBinHeights();
171  stat[1] = 0;
172  if(equivalentBinEntries() != 0)
173  stat[1] = (sumBinHeights() * sumBinHeights()) / equivalentBinEntries();
174  stat[2] = m_sumwx;
175  double meanX = 0;
176  if(sumBinHeights() != 0) meanX = m_sumwx/ sumBinHeights();
177  stat[3] = (meanX*meanX + rmsX*rmsX) * sumBinHeights();
178  stat[4] = m_sumwy;
179  double meanY = 0;
180  if(sumBinHeights() != 0) meanY = m_sumwy/ sumBinHeights();
181  stat[5] = (meanY*meanY + rmsY*rmsY) * sumBinHeights();
182  stat[6] = 0;
183  m_rep->PutStats(&stat.front());
184  return true;
185 }
186 
187 void Gaudi::Histogram2D::copyFromAida(const IHistogram2D& h) {
188  // implement here the copy
189  delete m_rep;
190  const char* tit = h.title().c_str();
191  if (h.xAxis().isFixedBinning() && h.yAxis().isFixedBinning() )
192  m_rep = new TH2D(tit,tit,
193  h.xAxis().bins(),h.xAxis().lowerEdge(),h.xAxis().upperEdge(),
194  h.yAxis().bins(),h.yAxis().lowerEdge(),h.yAxis().upperEdge() );
195  else {
196  Edges eX, eY;
197  for (int i =0; i < h.xAxis().bins(); ++i)
198  eX.push_back(h.xAxis().binLowerEdge(i));
199  // add also upperedges at the end
200  eX.push_back(h.xAxis().upperEdge() );
201  for (int i =0; i < h.yAxis().bins(); ++i)
202  eY.push_back(h.yAxis().binLowerEdge(i));
203  // add also upperedges at the end
204  eY.push_back(h.yAxis().upperEdge() );
205  m_rep = new TH2D(tit,tit,eX.size()-1,&eX.front(),eY.size()-1,&eY.front());
206  }
207  m_xAxis.initialize(m_rep->GetXaxis(),true);
208  m_yAxis.initialize(m_rep->GetYaxis(),true);
209  m_rep->Sumw2();
210  m_sumEntries = 0;
211  m_sumwx = 0;
212  m_sumwy = 0;
213  // statistics
214  double sumw = h.sumBinHeights();
215  double sumw2 = 0;
216  if (h.equivalentBinEntries() != 0)
217  sumw2 = ( sumw * sumw ) /h.equivalentBinEntries();
218  double sumwx = h.meanX()*h.sumBinHeights();
219  double sumwx2 = (h.meanX()*h.meanX() + h.rmsX()*h.rmsX() )*h.sumBinHeights();
220  double sumwy = h.meanY()*h.sumBinHeights();
221  double sumwy2 = (h.meanY()*h.meanY() + h.rmsY()*h.rmsY() )*h.sumBinHeights();
222  double sumwxy = 0;
223 
224  // copy the contents in (AIDA underflow/overflow are -2,-1)
225  for (int i=-2; i < xAxis().bins(); ++i) {
226  for (int j=-2; j < yAxis().bins(); ++j) {
227  // root binning starts from one !
228  m_rep->SetBinContent(rIndexX(i), rIndexY(j), h.binHeight(i,j) );
229  m_rep->SetBinError(rIndexX(i), rIndexY(j), h.binError(i,j) );
230  // calculate statistics
231  if ( i >= 0 && j >= 0) {
232  sumwxy += h.binHeight(i,j)*h.binMeanX(i,j)*h.binMeanY(i,j);
233  }
234  }
235  }
236  // need to do set entries after setting contents otherwise root will recalulate them
237  // taking into account how many time SetBinContents() has been called
238  m_rep->SetEntries(h.allEntries());
239  // fill stat vector
240  std::vector<double> stat(11);
241  stat[0] = sumw;
242  stat[1] = sumw2;
243  stat[2] = sumwx;
244  stat[3] = sumwx2;
245  stat[4] = sumwy;
246  stat[5] = sumwy2;
247  stat[6] = sumwxy;
248  m_rep->PutStats(&stat.front());
249 }
250 
251 typedef Gaudi::Histogram2D H2D;
253 
254 #ifdef __ICC
255 // re-enable icc remark #1572
256 #pragma warning(pop)
257 #endif

Generated at Wed Jan 30 2013 17:13:38 for Gaudi Framework, version v23r6 by Doxygen version 1.8.2 written by Dimitri van Heesch, © 1997-2004