ROOT  6.06/09
Reference Guide
SparseFit3.cxx
Go to the documentation of this file.
1 #include "TH2.h"
2 #include "TH3.h"
3 #include "TF2.h"
4 #include "TF3.h"
5 #include "TCanvas.h"
6 #include "TApplication.h"
7 #include "TMath.h"
8 
9 #include "Fit/SparseData.h"
10 #include "HFitInterface.h"
11 #include "Fit/Fitter.h"
12 #include "Math/WrappedMultiTF1.h"
13 
14 #include "TRandom.h"
15 
16 #include <iostream>
17 #include <iterator>
18 #include <algorithm>
19 
20 #include <list>
21 #include <vector>
22 
23 #include <cmath>
24 #include <cassert>
25 
26 using namespace std;
27 
28 bool showGraphics = false;
29 
30 double gaus2D(double *x, double *p)
31 {
32  return p[0]*TMath::Gaus(x[0],p[1],p[2]) * TMath::Gaus(x[1],p[3],p[4]);
33 }
34 
35 double gaus3D(double *x, double *p)
36 {
37  return p[0] * TMath::Gaus(x[0],p[1],p[2])
38  * TMath::Gaus(x[1],p[3],p[4])
39  * TMath::Gaus(x[2],p[5],p[6]);
40 }
41 
42 double pol2D(double *x, double *p)
43 {
44  return p[0]*x[0]+ p[1]*x[0]*x[0] + p[2]*x[1] + p[3]*x[1]*x[1] + p[4];
45 }
46 
47 ostream& operator << (ostream& out, ROOT::Fit::BinData& bd)
48 {
49  const unsigned int ndim( bd.NDim() );
50  const unsigned int npoints( bd.NPoints() );
51  for ( unsigned int i = 0; i < npoints; ++i )
52  {
53  double value, error;
54  const double *x = bd.GetPoint(i, value, error);
55  for ( unsigned int j = 0; j < ndim; ++j )
56  {
57  out << " x[" << j << "]: " << x[j];
58  }
59  out << " value: " << value;
60  out << " error: " << error;
61  out << endl;
62  }
63  return out;
64 }
65 
66 int findBin(ROOT::Fit::BinData& bd, const double *x)
67 {
68  const unsigned int ndim = bd.NDim();
69  const unsigned int npoints = bd.NPoints();
70 
71  for ( unsigned int i = 0; i < npoints; ++i )
72  {
73  double value1 = 0, error1 = 0;
74  const double *x1 = bd.GetPoint(i, value1, error1);
75  bool thisIsIt = true;
76  for ( unsigned int j = 0; j < ndim; ++j )
77  {
78  thisIsIt &= fabs(x1[j] - x[j]) < 1E-15;
79  }
80  if ( thisIsIt ) return i;
81  }
82 
83  cout << "NO ENCONTRADO!";
84  copy(x, x+ndim, ostream_iterator<double>(cout, " " ));
85  cout << endl;
86 
87  return -1;
88 }
89 
91 {
92  const unsigned int ndim = bd1.NDim();
93  const unsigned int npoints = bd1.NPoints();
94 
95  bool equals = true;
96 
97  cout << "Equals" << endl;
98 
99  for ( unsigned int i = 0; i < npoints && equals; ++i )
100  {
101  double value1 = 0, error1 = 0;
102  const double *x1 = bd1.GetPoint(i, value1, error1);
103 
104  int bin = findBin(bd2, x1);
105 
106  double value2 = 0, error2 = 0;
107  const double *x2 = bd2.GetPoint(bin, value2, error2);
108 
109  equals &= ( value1 == value2 );
110  cout << " v: " << equals;
111  equals &= ( error1 == error2 );
112  cout << " e: " << equals;
113  for ( unsigned int j = 0; j < ndim; ++j )
114  {
115  equals &= fabs(x1[j] - x2[j]) < 1E-15;
116  cout << " x[" << j << "]: " << equals;
117  }
118 
119  cout << " bd1: ";
120  std::copy(x1, &x1[ndim], ostream_iterator<double>(cout, " "));
121  cout << " value:" << value1;
122  cout << " error:" << error1;
123 
124  cout << " bd2: ";
125  std::copy(x2, &x2[ndim], ostream_iterator<double>(cout, " "));
126  cout << " value:" << value2;
127  cout << " error:" << error2;
128 
129  cout << " equals: " << equals;
130 
131  cout << endl;
132  }
133 
134  return equals;
135 }
136 
137 void fit3DHist()
138 {
139  vector<double> min(3); min[0] = 0.; min[1] = 0.; min[2] = 0.;
140  vector<double> max(3); max[0] = 10.; max[1] = 10.; max[2] = 10.;
141  vector<int> nbins(3); nbins[0] = 10; nbins[1] = 10; nbins[2] = 10;
142 
143  TH3D* h1 = new TH3D("3D Original Hist Fit", "h1-title",
144  nbins[0], min[0], max[0],
145  nbins[1], min[1], max[1],
146  nbins[2], min[2], max[2]);
147  TH3D* h2 = new TH3D("3D Blanked Hist Fit", "h1-title",
148  nbins[0], min[0], max[0],
149  nbins[1], min[1], max[1],
150  nbins[2], min[2], max[2]);
151 
152  TF3* f1 = new TF3("func3D", gaus3D,
153  min[0],max[0],
154  min[1],max[1],
155  min[2],max[2],
156  7);
157  double initialPars[] = {20,5,2,5,1,5,2};
158 // TF1* f1 = new TF1("func3D", pol2D,
159 // min[0],max[0], min[1], max[1], 5);
160 // double initialPars[] = {1,0.,0.5,0.,0.};
161  f1->SetParameters(initialPars);
162 // f1->FixParameter(1,0.);
163 // f1->FixParameter(3,0.);
164 
165 
166  for (int ix=1; ix <= h1->GetNbinsX(); ++ix) {
167  for (int iy=1; iy <= h1->GetNbinsY(); ++iy) {
168  for (int iz=1; iz <= h1->GetNbinsZ(); ++iz) {
169  double x = h1->GetXaxis()->GetBinCenter(ix);
170  double y = h1->GetYaxis()->GetBinCenter(iy);
171  double z = h1->GetZaxis()->GetBinCenter(iz);
172 
173  h1->SetBinContent( ix, iy, iz, gRandom->Poisson( f1->Eval(x,y,z) ) );
174  }
175  }
176  }
177 
178 ///////////////// CREATE THE SPARSE DATA
179  cout << "Retrieving the Sparse Data Structure" << endl;
180  ROOT::Fit::SparseData d(min,max);
181  ROOT::Fit::FillData(d, h1, 0);
183  d.GetBinData(bd);
184 
185 
186  cout << "Filling second histogram" << endl;
187  for ( unsigned int i = 0; i < bd.NPoints(); ++i )
188  {
189  const double* x;
190  double value, error;
191  x = bd.GetPoint(i, value, error);
192  value = (value)?value:-10;
193  h2->Fill(x[0], x[1], x[2], value);
194  }
195 
196 
197  ///////////////// FITS
198  // Fit preparation
199  bool ret;
200  ROOT::Fit::Fitter fitter;
203  fitter.Config().SetMinimizer("Minuit2");
204  //fitter.Config().MinimizerOptions().SetPrintLevel(3);
205 
206  /////////////////////////////////////////////////////////////////////////
207  cout << "\n ******* Fit with Original BinData *******" << endl;
208  ROOT::Fit::BinData bdOriginal;
209  ROOT::Fit::FillData(bdOriginal, h1, 0);
210  ret = fitter.LikelihoodFit(bdOriginal, if1, true);
211  fitter.Result().Print(std::cout);
212  if (!ret)
213  std::cout << "Fit Failed " << std::endl;
214 
215  /////////////////////////////////////////////////////////////////////////
216  cout << "\n ******* Fit with Original BinData with Zeros (integrate function in bins)*******" << endl;
218  opt.fUseEmpty = true;
219  opt.fIntegral = true;
220  ROOT::Fit::BinData bdOriginalWithCeros(opt);
221  ROOT::Fit::FillData(bdOriginalWithCeros, h1, 0);
222  ret = fitter.LikelihoodFit(bdOriginalWithCeros, if1, true);
223  fitter.Result().Print(std::cout);
224  if (!ret)
225  std::cout << "Fit Failed " << std::endl;
226 
227  /////////////////////////////////////////////////////////////////////////
228  cout << "\n ******* Fit with BinData and NoZeros *******" << endl;
229  ROOT::Fit::BinData bdNoCeros;
230  d.GetBinDataNoZeros(bdNoCeros);
231  ret = fitter.LikelihoodFit(bdNoCeros, if1, true);
232  fitter.Result().Print(std::cout);
233  if (!ret)
234  std::cout << "Fit Failed " << std::endl;
235 
236 // cout << "bdOriginal:\n" << bdOriginal << endl;
237 // cout << "bdNoCeros:\n" << *bdNoCeros << endl;
238 // cout << "Equals: " << (bdOriginal == *bdNoCeros) << endl;
239  /////////////////////////////////////////////////////////////////////////
240  cout << "\n ******* Fit with BinData with Zeros ******* (integrate function in bins)" << endl;
241  ROOT::Fit::BinData bdWithCeros(opt);
242  d.GetBinDataIntegral(bdWithCeros);
243  ret = fitter.LikelihoodFit(bdWithCeros, if1, true);
244  fitter.Result().Print(std::cout);
245  if (!ret)
246  std::cout << "Fit Failed " << std::endl;
247 
248 // cout << "bdOriginal:\n" << bdOriginal << endl;
249 // cout << "bdWithCeros:\n" << bdWithCeros << endl;
250 // cout << "Equals: " << (bdOriginal == bdWithCeros) << endl;
251 
252  /////////////////////////////////////////////////////////////////////////
253 
254  if (showGraphics) {
255 
256  TCanvas* c = new TCanvas("Histogram 3D");
257  c->Divide(1,2);
258  c->cd(1);
259  h1->Draw("lego");
260  c->cd(2);
261  h2->Draw("lego");
262  c->Update();
263  }
264 }
265 
266 void fit2DHist()
267 {
268  vector<double> min(2); min[0] = 0.; min[1] = 0.;
269  vector<double> max(2); max[0] = 10.; max[1] = 10.;
270  vector<int> nbins(2); nbins[0] = 10; nbins[1] = 10;
271 
272  TH2D* h1 = new TH2D("2D Original Hist Fit", "h1-title", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
273  TH2D* h2 = new TH2D("2D Blanked Hist Fit", "h1-title", nbins[0], min[0], max[0], nbins[1], min[1], max[1]);
274 
275  TF2* f2 = new TF2("func2D", gaus2D,
276  min[0],max[0], min[1], max[1], 5);
277  double initialPars[] = {20,5,2,5,1};
278 // TF2* f2 = new TF2("func2D", pol2D,
279 // min[0],max[0], min[1], max[1], 5);
280 // double initialPars[] = {1,0.,0.5,0.,0.};
281  f2->SetParameters(initialPars);
282 // f2->FixParameter(1,0.);
283 // f2->FixParameter(3,0.);
284 
285 
286  for (int ix=1; ix <= h1->GetNbinsX(); ++ix) {
287  for (int iy=1; iy <= h1->GetNbinsY(); ++iy) {
288  double x= h1->GetXaxis()->GetBinCenter(ix);
289  double y= h1->GetYaxis()->GetBinCenter(iy);
290 
291  h1->SetBinContent( ix, iy, gRandom->Poisson( f2->Eval(x,y) ) );
292  }
293  }
294 
295 ///////////////// CREATE THE SPARSE DATA
296  cout << "Retrieving the Sparse Data Structure" << endl;
297  ROOT::Fit::SparseData d(min,max);
298  ROOT::Fit::FillData(d, h1, 0);
300  d.GetBinData(bd);
301 
302 
303  cout << "Filling second histogram" << endl;
304  for ( unsigned int i = 0; i < bd.NPoints(); ++i )
305  {
306  const double* x;
307  double value, error;
308  x = bd.GetPoint(i, value, error);
309  value = (value)?value:-10;
310  h2->Fill(x[0], x[1], value);
311  }
312 
313  ///////////////// FITS
314  // Fit preparation
315  bool ret;
316  ROOT::Fit::Fitter fitter;
319  fitter.Config().SetMinimizer("Minuit2");
320  //fitter.Config().MinimizerOptions().SetPrintLevel(3);
321 
322  /////////////////////////////////////////////////////////////////////////
323  cout << "\n ******* Fit with Original BinData *******" << endl;
324  ROOT::Fit::BinData bdOriginal;
325  ROOT::Fit::FillData(bdOriginal, h1, 0);
326  ret = fitter.LikelihoodFit(bdOriginal, if2, true);
327  fitter.Result().Print(std::cout);
328  if (!ret)
329  std::cout << "Fit Failed " << std::endl;
330 
331  /////////////////////////////////////////////////////////////////////////
332  cout << "\n ******* Fit with Original BinData with Zeros*******" << endl;
334  opt.fUseEmpty = true;
335  opt.fIntegral = true;
336  ROOT::Fit::BinData bdOriginalWithCeros(opt);
337  ROOT::Fit::FillData(bdOriginalWithCeros, h1, 0);
338  ret = fitter.LikelihoodFit(bdOriginalWithCeros, if2, true);
339  fitter.Result().Print(std::cout);
340  if (!ret)
341  std::cout << "Fit Failed " << std::endl;
342 
343  /////////////////////////////////////////////////////////////////////////
344  cout << "\n ******* Fit with BinData and NoZeros *******" << endl;
345  ROOT::Fit::BinData bdNoCeros;
346  d.GetBinDataNoZeros(bdNoCeros);
347  ret = fitter.LikelihoodFit(bdNoCeros, if2, true);
348  fitter.Result().Print(std::cout);
349  if (!ret)
350  std::cout << "Fit Failed " << std::endl;
351 
352 // cout << "bdOriginal:\n" << bdOriginal << endl;
353 // cout << "bdNoCeros:\n" << *bdNoCeros << endl;
354 // cout << "Equals: " << (bdOriginal == *bdNoCeros) << endl;
355  /////////////////////////////////////////////////////////////////////////
356  cout << "\n ******* Fit with BinData with Zeros *******" << endl;
357  ROOT::Fit::BinData bdWithCeros(opt);
358  d.GetBinDataIntegral(bdWithCeros);
359  ret = fitter.LikelihoodFit(bdWithCeros, if2, true);
360  fitter.Result().Print(std::cout);
361  if (!ret)
362  std::cout << "Fit Failed " << std::endl;
363 
364 // cout << "bdOriginal:\n" << bdOriginal << endl;
365 // cout << "bdWithCeros:\n" << *bdWithCeros << endl;
366 // cout << "Equals: " << (bdOriginal == *bdWithCeros) << endl;
367 
368  /////////////////////////////////////////////////////////////////////////
369 
370  if (showGraphics) {
371  TCanvas* c = new TCanvas("Histogram 2D");
372  c->Divide(2,2);
373  c->cd(1);
374  h1->Draw("colz");
375  c->cd(2);
376  h1->Draw("text");
377  c->cd(3);
378  h2->Draw("colz");
379  c->cd(4);
380  h2->Draw("text");
381  c->Update();
382  }
383 }
384 
385 // void fit1DHist()
386 // {
387 // vector<double> min(1);
388 // min[0] = 0.;
389 
390 // vector<double> max(1);
391 // max[0] = 10.;
392 
393 // vector<int> nbins(1);
394 // nbins[0] = 10;
395 
396 // TH1D* h1 = new TH1D("1D Original Hist Fit", "h1-Original", nbins[0], min[0], max[0]);
397 // TH1D* h2 = new TH1D("1D Blanked Hist Fit", "h1-Blanked", nbins[0], min[0], max[0]);
398 
399 // TF1* f1 = new TF1("MyGaus", "[0]*TMath::Gaus([1],[2])", min[0], max[0]);
400 // f1->SetParameters(10., 5., 2.);
401 
402 // h1->FillRandom("MyGaus",1000);
403 
404 // cout << "Retrieving the Sparse Data Structure" << endl;
405 // ROOT::Fit::SparseData d(h1);
406 // ROOT::Fit::FillData(d, h1, 0);
407 // ROOT::Fit::BinData* bd = d.GetBinData();
408 
409 // cout << "Filling second histogram" << endl;
410 // for ( unsigned int i = 0; i < bd->NPoints(); ++i)
411 // {
412 // const double* x;
413 // double value, error;
414 // x = bd->GetPoint(i, value, error);
415 // value = (value)?value:-10;
416 // h2->Fill(x[0], value);
417 // }
418 
419 // TCanvas* c = new TCanvas("Histogram 2D");
420 // c->Divide(1,2);
421 // c->cd(1);
422 // h1->Draw("lego2Z");
423 // c->cd(2);
424 // h2->Draw("lego2Z");
425 
426 // // Fit preparation
427 // bool ret;
428 // ROOT::Fit::Fitter fitter;
429 // ROOT::Math::WrappedMultiTF1 wf1(*f1);
430 // fitter.Config().SetMinimizer("TMinuit");
431 
432 // cout << "\n ******* Chi2Fit with Original BinData *******" << endl;
433 // ROOT::Fit::BinData bdOriginal;
434 // ROOT::Fit::FillData(bdOriginal, h1, 0);
435 // ret = fitter.Fit(bdOriginal, wf1);
436 // fitter.Result().Print(std::cout);
437 // if (!ret)
438 // std::cout << "Fit Failed " << std::endl;
439 
440 // cout << "\n ******* Chi2Fit with BinData and NoCeros *******" << endl;
441 // ROOT::Fit::BinData* bdNoCeros = d.GetBinDataNoCeros();
442 
443 // cout << "bdOriginal:\n" << bdOriginal << endl;
444 // cout << "bdNoCeros:\n" << *bdNoCeros << endl;
445 
446 // cout << "Equals: " << (bdOriginal == *bdNoCeros) << endl;
447 
448 // ret = fitter.Fit(*bdNoCeros, wf1);
449 // fitter.Result().Print(std::cout);
450 // if (!ret)
451 // std::cout << "Fit Failed " << std::endl;
452 
453 
454 // delete bd;
455 // }
456 
457  int main(int argc, char** argv)
458 {
459 
460  bool do3dfit = false;
461  // Parse command line arguments
462  for (Int_t i=1 ; i<argc ; i++) {
463  std::string arg = argv[i] ;
464  if (arg == "-g") {
465  showGraphics = true;
466  }
467  if (arg == "-v") {
468  showGraphics = true;
469  //verbose = true;
470  }
471  if (arg == "-3d") {
472  do3dfit = true;
473  }
474  if (arg == "-h") {
475  cerr << "Usage: " << argv[0] << " [-g] [-v]\n";
476  cerr << " where:\n";
477  cerr << " -g : graphics mode\n";
478  cerr << " -v : verbose mode\n";
479  cerr << " -3d : 3d fit";
480  cerr << endl;
481  return -1;
482  }
483  }
484 
485  TApplication* theApp = 0;
486 
487  if ( showGraphics )
488  theApp = new TApplication("App",&argc,argv);
489 
490  fit2DHist();
491  if (do3dfit) fit3DHist();
492 // // fit1DHist();
493 
494  if ( showGraphics ) {
495  theApp->Run();
496  delete theApp;
497  theApp = 0;
498  }
499 
500  return 0;
501 }
502 
const FitConfig & Config() const
access to the fit configuration (const method)
Definition: Fitter.h:382
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:439
static Vc_ALWAYS_INLINE int_v min(const int_v &x, const int_v &y)
Definition: vector.h:433
THist< 2, double > TH2D
Definition: THist.h:320
Class to Wrap a ROOT Function class (like TF1) in a IParamMultiFunction interface of multi-dimensions...
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH3.cxx:3495
int main(int argc, char **argv)
Definition: SparseFit3.cxx:457
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
int Int_t
Definition: RtypesCore.h:41
int nbins[3]
int findBin(ROOT::Fit::BinData &bd, const double *x)
Definition: SparseFit3.cxx:66
int equals(Double_t n1, Double_t n2, double ERRORLIMIT=1.E-10)
Definition: UnitTesting.cxx:24
virtual Int_t GetNbinsX() const
Definition: TH1.h:296
STL namespace.
void fit2DHist()
Definition: SparseFit3.cxx:266
bool operator==(ROOT::Fit::BinData &bd1, ROOT::Fit::BinData &bd2)
Definition: SparseFit3.cxx:90
static const double x2[5]
Double_t x[n]
Definition: legend1.C:17
virtual void Run(Bool_t retrn=kFALSE)
Main application eventloop. Calls system dependent eventloop via gSystem.
double pol2D(double *x, double *p)
Definition: SparseFit3.cxx:42
const FitResult & Result() const
get fit result
Definition: Fitter.h:354
TH1F * h1
Definition: legend1.C:5
char * out
Definition: TBase64.cxx:29
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
void SetMinimizer(const char *type, const char *algo=0)
set minimizer type
Definition: FitConfig.h:152
void FillData(BinData &dv, const TH1 *hist, TF1 *func=0)
fill the data vector from a TH1.
void fit3DHist()
Definition: SparseFit3.cxx:137
void GetBinDataNoZeros(BinData &) const
Definition: SparseData.cxx:337
DataOptions : simple structure holding the options on how the data are filled.
Definition: DataOptions.h:28
IParamFunction interface (abstract class) describing multi-dimensional parameteric functions It is a ...
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
3-D histogram with a double per channel (see TH1 documentation)}
Definition: TH3.h:309
A 3-Dim function with parameters.
Definition: TF3.h:30
Fitter class, entry point for performing all type of fits.
Definition: Fitter.h:94
void GetBinDataIntegral(BinData &) const
Definition: SparseData.cxx:319
Double_t E()
Definition: TMath.h:54
virtual Int_t GetNbinsZ() const
Definition: TH1.h:298
TAxis * GetYaxis()
Definition: TH1.h:320
Class describing the binned data sets : vectors of x coordinates, y values and optionally error on y ...
Definition: BinData.h:61
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH3.cxx:275
A 2-Dim function with parameters.
Definition: TF2.h:33
R__EXTERN TRandom * gRandom
Definition: TRandom.h:62
bool LikelihoodFit(const BinData &data, bool extended=true)
Binned Likelihood fit.
Definition: Fitter.h:181
Double_t Gaus(Double_t x, Double_t mean=0, Double_t sigma=1, Bool_t norm=kFALSE)
Calculate a gaussian function with mean and sigma.
Definition: TMath.cxx:453
ostream & operator<<(ostream &out, ROOT::Fit::BinData &bd)
Definition: SparseFit3.cxx:47
THist< 3, double > TH3D
Definition: THist.h:326
The Canvas class.
Definition: TCanvas.h:48
unsigned int NPoints() const
return number of fit points
Definition: BinData.h:442
static const double x1[5]
void GetBinData(BinData &) const
Definition: SparseData.cxx:297
Double_t y[n]
Definition: legend1.C:17
void Print(std::ostream &os, bool covmat=false) const
print the result and optionaly covariance matrix and correlations
Definition: FitResult.cxx:446
static Vc_ALWAYS_INLINE int_v max(const int_v &x, const int_v &y)
Definition: vector.h:440
TAxis * GetZaxis()
Definition: TH1.h:321
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:449
virtual Int_t GetNbinsY() const
Definition: TH1.h:297
double gaus3D(double *x, double *p)
Definition: SparseFit3.cxx:35
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1077
double f2(const double *x)
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
Evaluate this function.
Definition: TF1.cxx:1185
TF1 * f1
Definition: legend1.C:11
This class creates the ROOT Application Environment that interfaces to the windowing system eventloop...
Definition: TApplication.h:45
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH2.cxx:2689
const double * GetPoint(unsigned int ipoint, double &value) const
retrieve at the same time a pointer to the coordinate data and the fit value More efficient than call...
Definition: BinData.h:304
virtual Int_t Poisson(Double_t mean)
Generates a random integer N according to a Poisson law.
Definition: TRandom.cxx:362
virtual void Update()
Update canvas pad buffers.
Definition: TCanvas.cxx:2179
Int_t Fill(Double_t)
Invalid Fill method.
Definition: TH2.cxx:285
bool showGraphics
Definition: SparseFit3.cxx:28
float value
Definition: math.cpp:443
double gaus2D(double *x, double *p)
Definition: SparseFit3.cxx:30
TAxis * GetXaxis()
Definition: TH1.h:319
unsigned int NDim() const
return coordinate data dimension
Definition: BinData.h:452
2-D histogram with a double per channel (see TH1 documentation)}
Definition: TH2.h:297