ATLAS Offline Software
Loading...
Searching...
No Matches
Resplot.h
Go to the documentation of this file.
1/* emacs: this is -*- c++ -*- */
21
22
23#ifndef RESPLOT_RESPLOT_H
24#define RESPLOT_RESPLOT_H
25
26#include <stdlib.h>
27
28#include <iostream>
29#include <string>
30#include <vector>
31
32#include "TROOT.h"
33#include "TH1D.h"
34#include "TH2D.h"
35#include "TF1.h"
36#include "TFile.h"
37#include "TString.h"
38
39
40#include "StatVal.h"
41
42
43// stupid root dictionary BOLLOX!! just because we inherit from a
44// TH2D, but I bet it *still* won't work because I use the
45// lovely "vector" class from the stl.
46#include "TObject.h"
47
48
49
50class Resplot : public TObject {
51
52private:
53
54 class ResException {
55 public:
56 ResException(const std::string& s) { std::cerr<<s<<std::endl; }
57 };
58
60
61public:
62
63 // NB, the destructor MUST NOT delete all the histograms, since root will
64 // when closing the file, and it will crash if we let the destructor
65 // try to clean up properly afterwards - LIKE WE SHOULD BE ABLE TO!
66 // NB: NB!! Actually, now we explicitly *remove* them from root managerment
67 // so we can do *proper* garbage collection, not like the
68 // alleged garbage collection in root
69
71 m_Set(false), m_Nentries(NULL),
72 m_mean(NULL), m_sigma(NULL), m_chi2(NULL),
73 m_h2d(NULL), m_h1d(NULL),
74 m_dir(NULL),
75 m_n_primary(0),
77 m_a_secondary(0.0),
78 m_b_secondary(0.0),
79 m_xaxis(""), m_yaxis(""), m_fitname(""), m_finalised(false),
80 m_uniform(true)
81 { }
82
83
84 Resplot(const std::string& name,
85 int n1, double a1, double b1,
86 int n2, double a2, double b2, const std::string& xaxis="") :
87 // TH2D( "2d", "2d", n1, a1, b1, n2, a2, b2),
88 m_Set(false), m_name(name),
89 m_Nentries(NULL),
90 m_mean(NULL), m_sigma(NULL), m_chi2(NULL),
91 m_h2d(NULL), m_h1d(NULL),
92 m_dir(NULL),
93 m_xaxis(xaxis), m_yaxis(""), m_fitname(""), m_finalised(false),
94 m_uniform(true)
95 { Initialise(name, n1, a1, b1, n2, a2, b2); }
96
97
98 Resplot(const std::string& name,
99 int n1, const double* a1,
100 int n2, double a2, double b2, const std::string& xaxis="") :
101 // TH2D( "2d", "2d", n1, a1, n2, a2, b2),
102 m_Set(false), m_name(name),
103 m_Nentries(NULL),
104 m_mean(NULL), m_sigma(NULL), m_chi2(NULL),
105 m_h2d(NULL), m_h1d(NULL),
106 m_dir(NULL),
107 m_xaxis(xaxis), m_yaxis(""), m_fitname(""), m_finalised(false),
108 m_uniform(true)
109 { Initialise(name, n1, a1, n2, a2, b2); }
110
111 Resplot(const std::string& name,
112 const std::vector<double>& a,
113 int n2, double a2, double b2, const std::string& xaxis="") :
114 m_Set(false), m_name(name),
115 m_Nentries(NULL),
116 m_mean(NULL), m_sigma(NULL), m_chi2(NULL),
117 m_h2d(NULL), m_h1d(NULL),
118 m_dir(NULL),
119 m_xaxis(xaxis), m_yaxis(""), m_slices(), m_fitname(""), m_finalised(false),
120 m_uniform(true)
121 { Initialise(name, a, n2, a2, b2); }
122
123
124 Resplot(const std::string& name,
125 const std::vector<double>& a,
126 const std::vector<double>& b,
127 const std::string& xaxis="") :
128 m_Set(false), m_name(name),
129 m_Nentries(NULL),
130 m_mean(NULL), m_sigma(NULL), m_chi2(NULL),
131 m_h2d(NULL), m_h1d(NULL),
132 m_dir(NULL),
133 m_xaxis(xaxis), m_yaxis(""), m_slices(), m_fitname(""), m_finalised(false),
134 m_uniform(false)
135 { Initialise(name, a, b); }
136
137
138
139 Resplot(const Resplot& r, const std::string& tag="") :
140 // TH2D( (TH2D)r ),
141 TObject(r),
142 m_Set(r.m_Set), m_name(r.m_name+tag),
143 m_Nentries(new TH1D(*r.m_Nentries)),
144 m_mean(new TH1D(*r.m_mean)), m_sigma(new TH1D(*r.m_sigma)), m_chi2(new TH1D(*r.m_chi2)),
145 m_h2d(new TH2D(*r.m_h2d)), m_h1d(new TH1D(*r.m_h1d)),
146 m_dir(NULL),
147 m_n_primary(r.m_n_primary), // a_primary(r.a_primary), b_primary(r.b_primary),
150 m_uniform(true)
151 {
152 // TH2D::SetDirectory(0);
153 skip(m_h2d);
154 skip(m_h1d);
155 skip(m_mean);
156 skip(m_sigma);
157 skip(m_chi2);
159 }
160
161 Resplot(const std::string& name, TH2* hin, const bool flip=false) :
162 m_Set(false), m_name(name),
163 m_Nentries(NULL),
164 m_mean(NULL), m_sigma(NULL), m_chi2(NULL),
165 m_h2d(NULL), m_h1d(NULL),
166 m_dir(NULL),
167 m_xaxis(""), m_yaxis(""), m_fitname(""), m_finalised(false),
168 m_uniform(true)
169 {
170
171 TH2D* h = 0;
172
173 if ( flip ) h = rotate( hin );
174 else h = static_cast<TH2D*>(hin);
175
176 // get the major bin edges
177 std::vector<double> a1;
178
179 TAxis* axe = h->GetXaxis();
180
181 for ( int i=1 ; i<=axe->GetNbins()+1 ; i++ ) a1.push_back(axe->GetBinLowEdge(i));
182
183 // get the y bin limits
184 TH1D* s = h->ProjectionY("_duff", 1, 1, "e");
185
186 int n2 = s->GetNbinsX();
187 double a2 = s->GetBinLowEdge(1);
188 double b2 = s->GetBinLowEdge(s->GetNbinsX()+1);
189
190 delete s;
191
192 Initialise(name, a1, n2, a2, b2 );
193
194 // replace default histograms
195 delete m_h2d;
196 m_h2d = (TH2D*)h->Clone("2d");
197 skip(m_h2d);
198
199 delete m_h1d;
200 m_h1d = h->ProjectionY("1d", 1, h->GetNbinsX(), "e");
201 m_h1d->SetTitle("1d");
202 skip(m_h1d);
203 }
204
205
206 // hmmmm, this constructor doesn't work for some
207 // reason. Can't work out why, must be due to kak root directory
208 // object ownership, I shouldn't wonder
209 Resplot(const std::string& name) :
210 m_name(name)
211 // m_Nentries ( new TH1D(*(TH1D*)gDirectory->Get((name+"/fractional uncertainty").c_str())) ),
212 // m_mean ( new TH1D(*(TH1D*)gDirectory->Get((name+"/1d").c_str())) ),
213 // m_sigma ( new TH1D(*(TH1D*)gDirectory->Get((name+"/sigma").c_str())) ),
214 // m_chi2 ( new TH1D(*(TH1D*)gDirectory->Get((name+"/chi2").c_str())) ),
215 // m_h2d ( new TH2D(*(TH2D*)gDirectory->Get((name+"/2d").c_str())) ),
216 // m_h1d ( new TH1D(*(TH1D*)gDirectory->Get((name+"/1d").c_str())) )
217 {
218
219 m_Nentries = new TH1D(*(TH1D*)gDirectory->Get((name+"/fractional uncertainty").c_str()));
220 m_mean = new TH1D(*(TH1D*)gDirectory->Get((name+"/mean").c_str()));
221 m_sigma = new TH1D(*(TH1D*)gDirectory->Get((name+"/sigma").c_str()));
222 // m_sigma = new TH1D(*(TH1D*)gDirectory->Get((name+"/RMS").c_str()));
223 m_chi2 = new TH1D(*(TH1D*)gDirectory->Get((name+"/chi2").c_str()));
224 m_h2d = new TH2D(*(TH2D*)gDirectory->Get((name+"/2d").c_str()));
225 m_h1d = new TH1D(*(TH1D*)gDirectory->Get((name+"/1d").c_str()));
226
227 skip(m_h2d);
228 skip(m_h1d);
229 skip(m_mean);
230 skip(m_sigma);
231 skip(m_chi2);
233
234 m_n_primary = m_h2d->GetNbinsX();
235
236 m_finalised = true;
237
238 m_uniform = true;
239
240 // std::cout << "Resplot(std::string)" << name << "\t m_n_primary " << m_n_primary << std::endl;
241
242 // deep copy
243 for ( int i=1 ; i<=m_h2d->GetNbinsX() ; i++ ) {
244 char slicename[512];
245 sprintf(slicename, "/slices/slice[%03d]", i);
246 // std::cout << i << " " << (name+slicename) << std::endl;
247 m_slices.push_back((TH1D*)gDirectory->Get((name+slicename).c_str()));
248 skip(m_slices.back()); //->SetDirectory(0);
249 }
250
251 }
252
253
254 virtual ~Resplot();
255
256
257 void Initialise(const std::string& name,
258 int n1, double a1, double b1,
259 int n2, double a2, double b2);
260
261 void Initialise(const std::string& name,
262 int n1, const double* a1,
263 int n2, double a2, double b2);
264
265 void Initialise(const std::string& name,
266 int n1, const double* a1,
267 int n2, const double* a2);
268
269 void Initialise(const std::string& name,
270 const std::vector<double>& a,
271 int n2, double a2, double b2);
272
273 void Initialise(const std::string& name,
274 const std::vector<double>& a,
275 const std::vector<double>& b );
276
277
278 // Fill the helper histograms
279 int Fill(double x, double y, double w=1) {
280 // m_h2d->Fill(x, y, w);
281 // TH2D::Fill(x, y, w);
282 m_h2d->Fill(x, y, w);
283 m_h1d->Fill(y, w);
284 return 0;
285 }
286
287
288 // Get the efficiency for a single slice of the major ordinate
289 StatVal GetEfficiency(TH1D* h, double lower, double upper);
290
291
292 // Get the global efficiency
293 StatVal GetGlobalEfficiency(double lower, double upper) {
294 return GetEfficiency(m_h1d, lower, upper);
295 }
296
297
298 // Get all efficiencies for all the major bins
299 // within the range lower .. upper
300 TH1D* GetEfficiencies(const std::string& hname, double lower, double upper);
301
302 // Get all efficiencies for all the major bins
303 // within the range += Nsigma*sigma about the fitted
304 // mean, where the mean and sigma are taken from the fit,
305 // either gaussian of cauchy
306 TH1D* GetEfficiencies(const std::string& hname, double Nsigma);
307
308 // set the name of the x axis variable
309 void SetXAxis(const std::string& xaxis) {
310 m_xaxis = xaxis;
311 m_mean->SetXTitle(m_xaxis.c_str());
312 m_sigma->SetXTitle(m_xaxis.c_str());
313 m_chi2->SetXTitle(m_xaxis.c_str());
314 m_Nentries->SetXTitle(m_xaxis.c_str());
315 m_h2d->SetXTitle(m_xaxis.c_str());
316 m_h1d->SetXTitle(m_xaxis.c_str());
317 }
318
319 // and the y axis variable
320 void SetYAxis(const std::string& yaxis) {
321 m_yaxis = yaxis;
322 m_h2d->SetXTitle(m_yaxis.c_str());
323 m_h1d->SetXTitle(m_yaxis.c_str());
324 }
325
326 // get the limits of the slices in the x coordinate
327 std::vector<double> GetPrimaryLimits() const {
328 std::vector<double> d;
329 for ( int i=1 ; i<=m_mean->GetNbinsX()+1 ; i++ ) d.push_back(m_mean->GetBinLowEdge(i));
330 return d;
331 }
332
333 // get the central values for each slice in the x coordinate
334 std::vector<double> GetPrimaryCentres() const {
335 std::vector<double> d;
336 for ( int i=1 ; i<=m_mean->GetNbinsX() ; i++ ) d.push_back(m_mean->GetBinCenter(i));
337 return d;
338 }
339
340 const std::string& Name() const { return m_name; }
341
342 const std::string& Name(const std::string& s) { return m_name=s; }
343
344 // get the histogram of the mean of the distributions
345 TH1D* Mean() { return m_mean; }
346 const TH1D* Mean() const { return m_mean; }
347
348 // and the "width"
349 TH1D* Sigma() { return m_sigma; }
350
351 // and the chi2
352 TH1D* Chi2() { return m_chi2; }
353 const TH1D* Chi2() const { return m_chi2; }
354
355 // and the overall projection in y
356 TH1D* H1D() { return m_h1d; }
357 const TH1D* H1D() const { return m_h1d; }
358
359 // and the overall 2D histo
360 TH2D* H2D() { return m_h2d; }
361 const TH2D* H2D() const { return m_h2d; }
362
363 // return 1/sqrt(Nentries) approximate stat precision
364 TH1D* Uncertainty() { return m_Nentries; }
365 const TH1D* Uncertainty() const { return m_Nentries; }
366
367 bool finalised() const { return m_finalised; }
368
369 // and the fitted slices
370 const std::vector<TH1D*>& Slices() const { return m_slices; }
371
372 // and the name of the fitted function
373 // NB. this is useful for setting plot paramaters
374 // of the fit for plotting, such as the line colour etc
375 const std::string& FitName() const { return m_fitname; }
376
377 // Static fit functions!!!!
378 // NB. *Any* fit function can be used to fit the distributions
379 // as long as it has this form of method, it does not need
380 // to be a member of the class Resplot, but using static
381 // member functions here, gaurantees there are a few simple
382 // examples that will be safely encapsulated inside the
383 // class so there will be no interference with other functions
384
385 // fit Generic to histogram
386 static TF1* FitInternal(TH1D* s, double a=-999, double b=-999, TF1* f1=0);
387
388 // fit Gaussian to histogram
389 static TF1* FitGaussian(TH1D* s, double a=-999, double b=-999);
390
391 // fit Gaussian with radial dependence to histogram
392 static TF1* FitRGaussian(TH1D* s, double a=-999, double b=-999);
393
394 // or fit a Breit Wigner
395 static TF1* FitBreit(TH1D* s, double a=-999, double b=-999);
396
397 // or fit a Breit Wigner with radial dependence
398 static TF1* FitRBreit(TH1D* s, double a=-999, double b=-999);
399
400 // or fit an inverse tan turn on curve
401 static TF1* FitATan(TH1D* s, double a=-999, double b=-999);
402
403 // use histo mean and rms
404 static TF1* FitNull(TH1D* s, double a=-999, double b=-999);
405
406 // fit a poisson distribution
407 static TF1* FitPoisson(TH1D* s, double a=-999, double b=-999);
408
409 // fit a ax^b exp(c*x)+d distribution
410 static TF1* FitXExp(TH1D* s, double a=-999, double b=-999);
411
412 static TF1* FitLandauGauss(TH1D* s, double a=-999, double b=-999);
413
414 // fit Central Gaussian to histogram
415 static TF1* FitCentralGaussian(TH1D* s, double a=-999, double b=-999);
416
417 // fit histo mean and rms only to central 95% of histogram
418 static TF1* FitNull95(TH1D* s, double a=0, double b=0);
419 // static TF1* FitNull95(TH1D* s, double frac=0.95, bool useroot=false);
420 static TF1* FitNull95New(TH1D* s, double frac=0.95, bool useroot=false);
421 static TF1* FitNull95Obsolete(TH1D* s, double frac=0.95, bool useroot=false);
422
423 static TF1* FitNull95Central(TH1D* s);
424
425
426
427 // Get the offsets and resolutions overall and for each slice
428 // int Finalise(double a=-999, double b=-999, TF1* (*func)(TH1D* s, double a, double b)=Resplot::FitGaussian);
429 int Finalise(double a=-999, double b=-999, TF1* (*func)(TH1D* s, double a, double b)=Resplot::FitGaussian);
430
431 int Finalise(TF1* (*func)(TH1D* s, double a, double b), double a=-999, double b=-999 ) {
432 return Finalise( a, b, func);
433 }
434
435 int Fit(TF1* (*func)(TH1D* s, double a, double b) ) {
436 if ( m_finalised ) return -1;
437 return Finalise( func );
438 }
439
440 int Refit(TF1* (*func)(TH1D* s, double a, double b) ) {
441 clear();
442 return Finalise( func );
443 }
444
445
446 TDirectory* Dir() { return m_dir; }
447
452 Int_t Write(const char* =0, Int_t =0, Int_t =0) const { return DWrite(); }
453 Int_t Write(const char* =0, Int_t =0, Int_t =0) { return DWrite(); }
454 Int_t DWrite(TDirectory* g=0) const;
455
457 static void AddDirectory(bool t=false) { s_mAddDirectoryStatus=t; }
458
459 void SetDirectory(TDirectory* =0) { }
460
461 static void setInterpolateFlag(bool b) { s_interpolate_flag = b; }
462
463
465
467
468 clear();
469
470 // std::cout << "Resplot::operator+= " << Name() << std::endl;
471
472 if ( r.H1D()->GetNbinsX()!=H1D()->GetNbinsX()) throw ResException("histogram limits do not match");
473 if ( r.H2D()->GetNbinsX()!=H2D()->GetNbinsX()) throw ResException("histogram limits do not match");
474 if ( r.H2D()->GetNbinsY()!=H2D()->GetNbinsY()) throw ResException("histogram limits do not match");
475
476 for ( int i=0 ; i<=r.H2D()->GetNbinsX()+1 ; i++ ) {
477 for ( int j=0 ; j<=r.H2D()->GetNbinsY()+1 ; j++ ) {
478 H2D()->SetBinContent( i, j, H2D()->GetBinContent(i,j) + r.H2D()->GetBinContent(i,j) );
479 H2D()->SetBinError( i, j, std::sqrt(H2D()->GetBinError(i,j)*H2D()->GetBinError(i,j) + r.H2D()->GetBinError(i,j)*r.H2D()->GetBinError(i,j) ) );
480 }
481 }
482
483 for ( int i=0 ; i<=r.H1D()->GetNbinsX()+1 ; i++ ) {
484 H1D()->SetBinContent( i, H1D()->GetBinContent(i) + r.H1D()->GetBinContent(i) );
485 H1D()->SetBinError( i, std::sqrt(H1D()->GetBinError(i)*H1D()->GetBinError(i) + r.H1D()->GetBinError(i)*r.H1D()->GetBinError(i) ) );
486 }
487
488 m_finalised = false;
489
490 // std::cout << "Resplot::operator+= done" << std::endl;
491
492 return *this;
493 }
494
495
496
497 Resplot& operator*=(double d) {
498
499 std::cout << "Resplot::operator*= " << Name() << " " << m_finalised << " " << m_slices.size() << std::endl;
500
501 clear();
502
503 std::cout << "Resplot::operator*= " << Name() << " " << m_finalised << " " << m_slices.size() << std::endl;
504 std::cout << "Resplot::operator*= " << Name() << "\tH2D " << H2D() << std::endl;
505
506 H2D()->DrawCopy();
507
508 std::cout << "Resplot::operator*= " << Name() << "\tDraw " << std::endl;
509
510 for ( int i=0 ; i<=H2D()->GetNbinsX()+1 ; i++ ) {
511 for ( int j=0 ; j<=H2D()->GetNbinsY()+1 ; j++ ) {
512 H2D()->SetBinContent( i, j, H2D()->GetBinContent(i,j)*d);
513 H2D()->SetBinError( i, j, H2D()->GetBinError(i,j)*d);
514 }
515 }
516
517
518 std::cout << "Resplot::operator*= " << Name() << "\tH1D " << H1D() << std::endl;
519
520 for ( int i=0 ; i<=H1D()->GetNbinsX()+1 ; i++ ) {
521 H1D()->SetBinContent( i, H1D()->GetBinContent(i)*d );
522 H1D()->SetBinError( i, H1D()->GetBinError(i)*d );
523 }
524
525 m_finalised = false;
526
527 std::cout << "Resplot::operator*= done" << std::endl;
528
529 return *this;
530 }
531
532 Resplot operator*(double d) const {
533 Resplot r(*this);
534 return (r*=d);
535 }
536
537 Resplot operator+(const Resplot& r) const {
538 Resplot r0(*this);
539 return (r0+=r);
540 }
541
542
543 void Integrate() {
544 if ( !m_finalised ) return;
545
548 }
549
550 void setUniform(bool t=false) { m_uniform=t; }
551
552 static const std::string& version() { return s_rversion; }
553 static bool setoldrms95(bool b) { return s_oldrms95=b; }
554 static bool setscalerms95(bool b) { return s_scalerms95=b; }
555 static bool setnofit(bool b) { return s_nofit=b; }
556
558 static TH2D* rotate(const TH2* h);
559
562 static TH2D* combine(const TH2* h, double x, int N);
564 static TH1D* combine(const TH1* h, double x, int N);
565
568 static TH2D* combine(const TH2* h, double epsilon);
569 static TH1D* combine(const TH1* h, double epsilon);
570
571 static std::vector<double> findbins(const TH1* h, double epsilon);
572 static TH1D* rebin(const TH1* h, const std::vector<double>& bins);
573
574 static std::vector<double> getbins(const TH1* h);
575
576private:
577
578 void MakeSlices() { }
579
580private:
581
582 // stop the root default histogram ownership nonsense
583 void skip(TH1D* t) { if ( t ) t->SetDirectory(0); }
584 void delskip(TH1D* t) { if ( t ) { t->SetDirectory(0); delete t; } }
585 void deletehist(TH1D* t) { if ( t ) delete t; }
586
587 void skip(TH2D* t) { if ( t ) t->SetDirectory(0); }
588 void delskip(TH2D* t) { if ( t ) { t->SetDirectory(0); delete t; } }
589 void deletehist(TH2D* t) { if ( t ) delete t; }
590
591
592 void SetPrimary(int n, const double* ) { m_n_primary=n; }
593 void SetPrimary(int n, double , double ) { m_n_primary=n; } // a_primary=a; b_primary=b; }
594 void SetSecondary(int n, double a, double b) { m_n_secondary=n; m_a_secondary=a; m_b_secondary=b; }
595
596
597 // stuff for calculating efficiencies
598
601
602
603 void AddEfficiency(TH1D* h, int i, double lower, double upper) {
604 StatVal v = GetGlobalEfficiency( lower, upper);
605 h->SetBinContent(i, v.value);
606 h->SetBinError(i, v.error);
607 }
608
609 void AddResolution(TH1D* h, int i) {
610 StatVal v = get_sigma();
611 h->SetBinContent(i, v.value);
612 h->SetBinError(i, v.error);
613 }
614
615 void AddMean(TH1D* h, int i) {
616 StatVal v = get_mean();
617 h->SetBinContent(i, v.value);
618 h->SetBinError(i, v.error);
619 }
620
621private:
622
623 void clear() {
624 // std::cout << "clear " << m_slices.size() << std::endl;
625 m_finalised = false;
626 for ( int i=m_slices.size() ; i-- ; ) if ( m_slices[i] ) delete m_slices[i];
627 m_slices.clear();
628 // std::cout << "cleared" << std::endl;
629 }
630
631private:
632
635
636 static const std::string s_rversion;
637
638 bool m_Set;
639
640 std::string m_name;
641
642 TH1D* m_Nentries = 0;
643 TH1D* m_mean = 0;
644 TH1D* m_sigma = 0;
645 TH1D* m_chi2 = 0;
646
647 TH2D* m_h2d = 0;
648 TH1D* m_h1d = 0;
649
652
653 TDirectory* m_dir = 0;
654
656 // double a_primary, b_primary;
657
660
661 std::string m_xaxis;
662 std::string m_yaxis;
663
664 std::vector<TH1D*> m_slices;
665 std::string m_fitname;
666
668 TDirectory* m_slicedir = 0;
669
671
674 static bool s_oldrms95;
675 static bool s_scalerms95;
676
677 static bool s_nofit;
678
679
681
682 ClassDef(Resplot, 1)
683
684};
685
686
687inline Resplot operator*(double d, Resplot r) { return (r*=d); }
688
689
690#endif /* RESPLOT_RESPLOT_H */
691
int upper(int c)
static Double_t a
static const std::vector< std::string > bins
Resplot operator*(double d, Resplot r)
Definition Resplot.h:687
#define y
#define x
Header file for AthHistogramAlgorithm.
ResException(const std::string &s)
Definition Resplot.h:56
static bool s_nofit
Definition Resplot.h:677
static TH2D * rotate(const TH2 *h)
flip the x and y axes
Definition Resplot.cxx:564
TH1D * H1D()
Definition Resplot.h:356
TH2D * m_h2d
Definition Resplot.h:647
std::vector< double > GetPrimaryLimits() const
Definition Resplot.h:327
const TH2D * H2D() const
Definition Resplot.h:361
static TF1 * FitRGaussian(TH1D *s, double a=-999, double b=-999)
Definition Resplot.cxx:1031
static std::vector< double > getbins(const TH1 *h)
Definition Resplot.cxx:921
std::string m_yaxis
Definition Resplot.h:662
const std::string & Name() const
Definition Resplot.h:340
static TF1 * FitCentralGaussian(TH1D *s, double a=-999, double b=-999)
Definition Resplot.cxx:1080
Int_t Write(const char *=0, Int_t=0, Int_t=0)
Definition Resplot.h:453
std::vector< TH1D * > m_slices
Definition Resplot.h:664
bool m_Set
Definition Resplot.h:638
static TF1 * FitInternal(TH1D *s, double a=-999, double b=-999, TF1 *f1=0)
Definition Resplot.cxx:1051
static TH1D * rebin(const TH1 *h, const std::vector< double > &bins)
Definition Resplot.cxx:891
static TF1 * FitNull95New(TH1D *s, double frac=0.95, bool useroot=false)
Definition Resplot.cxx:1687
void SetXAxis(const std::string &xaxis)
Definition Resplot.h:309
static const std::string s_rversion
Definition Resplot.h:636
static TF1 * FitGaussian(TH1D *s, double a=-999, double b=-999)
Definition Resplot.cxx:1025
static TF1 * FitRBreit(TH1D *s, double a=-999, double b=-999)
Definition Resplot.cxx:1036
int Finalise(double a=-999, double b=-999, TF1 *(*func)(TH1D *s, double a, double b)=Resplot::FitGaussian)
Definition Resplot.cxx:388
TDirectory * Dir()
Definition Resplot.h:446
std::string m_xaxis
Definition Resplot.h:661
const std::string & Name(const std::string &s)
Definition Resplot.h:342
Resplot & operator+=(const Resplot &r)
operators
Definition Resplot.h:466
static bool setscalerms95(bool b)
Definition Resplot.h:554
std::string m_name
Definition Resplot.h:640
StatVal m_g_mean
Definition Resplot.h:650
Int_t Write(const char *=0, Int_t=0, Int_t=0) const
Hooray, this stupidity is to overwride both the const and non-const TObject Write methods Fixme: shou...
Definition Resplot.h:452
Resplot(const std::string &name)
Definition Resplot.h:209
std::vector< double > GetPrimaryCentres() const
Definition Resplot.h:334
void AddMean(TH1D *h, int i)
Definition Resplot.h:615
TH1D * Chi2()
Definition Resplot.h:352
void SetYAxis(const std::string &yaxis)
Definition Resplot.h:320
static TF1 * FitBreit(TH1D *s, double a=-999, double b=-999)
Definition Resplot.cxx:1159
static bool AddDirectoryStatus()
Definition Resplot.h:456
Resplot(const std::string &name, int n1, double a1, double b1, int n2, double a2, double b2, const std::string &xaxis="")
Definition Resplot.h:84
void SetSecondary(int n, double a, double b)
Definition Resplot.h:594
TH1D * m_mean
Definition Resplot.h:643
void delskip(TH2D *t)
Definition Resplot.h:588
Resplot()
Definition Resplot.h:70
static TF1 * FitLandauGauss(TH1D *s, double a=-999, double b=-999)
Definition Resplot.cxx:2082
static TF1 * FitXExp(TH1D *s, double a=-999, double b=-999)
Definition Resplot.cxx:1966
Resplot(const std::string &name, const std::vector< double > &a, const std::vector< double > &b, const std::string &xaxis="")
Definition Resplot.h:124
int Refit(TF1 *(*func)(TH1D *s, double a, double b))
Definition Resplot.h:440
TH1D * m_h1d
Definition Resplot.h:648
TH1D * m_Nentries
Definition Resplot.h:642
static bool s_oldrms95
temporarily allow toggeling between old and new rms95 error estimates
Definition Resplot.h:674
static void AddDirectory(bool t=false)
Definition Resplot.h:457
void SetDirectory(TDirectory *=0)
Definition Resplot.h:459
TH2D * H2D()
Definition Resplot.h:360
Resplot(const std::string &name, TH2 *hin, const bool flip=false)
Definition Resplot.h:161
double m_b_secondary
Definition Resplot.h:659
StatVal m_g_sigma
Definition Resplot.h:651
const std::string & FitName() const
Definition Resplot.h:375
static const std::string & version()
Definition Resplot.h:552
static TF1 * FitNull(TH1D *s, double a=-999, double b=-999)
Definition Resplot.cxx:1288
void clear()
Definition Resplot.h:623
static TF1 * FitNull95(TH1D *s, double a=0, double b=0)
Definition Resplot.cxx:1673
static bool s_interpolate_flag
Definition Resplot.h:634
Resplot & operator*=(double d)
Definition Resplot.h:497
void AddEfficiency(TH1D *h, int i, double lower, double upper)
Definition Resplot.h:603
int m_n_primary
Definition Resplot.h:655
void SetPrimary(int n, double, double)
Definition Resplot.h:593
int Finalise(TF1 *(*func)(TH1D *s, double a, double b), double a=-999, double b=-999)
Definition Resplot.h:431
static TF1 * FitNull95Obsolete(TH1D *s, double frac=0.95, bool useroot=false)
Definition Resplot.cxx:1408
bool m_uniform
Definition Resplot.h:670
const TH1D * Uncertainty() const
Definition Resplot.h:365
static bool setoldrms95(bool b)
Definition Resplot.h:553
Resplot(const std::string &name, const std::vector< double > &a, int n2, double a2, double b2, const std::string &xaxis="")
Definition Resplot.h:111
TH1D * m_sigma
Definition Resplot.h:644
void skip(TH2D *t)
Definition Resplot.h:587
StatVal GetEfficiency(TH1D *h, double lower, double upper)
Definition Resplot.cxx:932
@ BINWIDTH
Definition Resplot.h:59
@ HISTOWIDTH
Definition Resplot.h:59
@ EMPTY
Definition Resplot.h:59
Resplot(const std::string &name, int n1, const double *a1, int n2, double a2, double b2, const std::string &xaxis="")
Definition Resplot.h:98
StatVal GetGlobalEfficiency(double lower, double upper)
Definition Resplot.h:293
void SetPrimary(int n, const double *)
Definition Resplot.h:592
void Initialise(const std::string &name, int n1, double a1, double b1, int n2, double a2, double b2)
Definition Resplot.cxx:194
Resplot(const Resplot &r, const std::string &tag="")
Definition Resplot.h:139
static TH2D * combine(const TH2 *h, double x, int N)
combine bins along the x axis after value x, combine each N bins into a single bin
Definition Resplot.cxx:606
const TH1D * Mean() const
Definition Resplot.h:346
void MakeSlices()
Definition Resplot.h:578
void skip(TH1D *t)
Definition Resplot.h:583
Resplot operator+(const Resplot &r) const
Definition Resplot.h:537
StatVal get_mean()
Definition Resplot.h:599
TH1D * m_chi2
Definition Resplot.h:645
void AddResolution(TH1D *h, int i)
Definition Resplot.h:609
static TF1 * FitPoisson(TH1D *s, double a=-999, double b=-999)
Definition Resplot.cxx:1926
double m_a_secondary
Definition Resplot.h:659
Int_t DWrite(TDirectory *g=0) const
boooo, need to use the stupid Int_t class because I foolishly gave this function the same name as TOb...
Definition Resplot.cxx:342
bool finalised() const
Definition Resplot.h:367
void deletehist(TH1D *t)
Definition Resplot.h:585
static ERROR s_ErrorSet
Definition Resplot.h:680
TH1D * Uncertainty()
Definition Resplot.h:364
const TH1D * H1D() const
Definition Resplot.h:357
static void setInterpolateFlag(bool b)
Definition Resplot.h:461
std::string m_fitname
Definition Resplot.h:665
int Fit(TF1 *(*func)(TH1D *s, double a, double b))
Definition Resplot.h:435
void Integrate()
Definition Resplot.h:543
Resplot operator*(double d) const
Definition Resplot.h:532
TDirectory * m_slicedir
Definition Resplot.h:668
TDirectory * m_dir
Definition Resplot.h:653
int Fill(double x, double y, double w=1)
Definition Resplot.h:279
static bool s_scalerms95
Definition Resplot.h:675
int m_n_secondary
Definition Resplot.h:658
static bool s_mAddDirectoryStatus
Definition Resplot.h:633
TH1D * Sigma()
Definition Resplot.h:349
const TH1D * Chi2() const
Definition Resplot.h:353
void setUniform(bool t=false)
Definition Resplot.h:550
void delskip(TH1D *t)
Definition Resplot.h:584
static TF1 * FitATan(TH1D *s, double a=-999, double b=-999)
Definition Resplot.cxx:1220
bool m_finalised
Definition Resplot.h:667
static TF1 * FitNull95Central(TH1D *s)
Definition Resplot.cxx:1754
StatVal get_sigma()
Definition Resplot.h:600
virtual ~Resplot()
Definition Resplot.cxx:174
TH1D * Mean()
Definition Resplot.h:345
TH1D * GetEfficiencies(const std::string &hname, double lower, double upper)
Definition Resplot.cxx:958
void deletehist(TH2D *t)
Definition Resplot.h:589
const std::vector< TH1D * > & Slices() const
Definition Resplot.h:370
static bool setnofit(bool b)
Definition Resplot.h:555
static std::vector< double > findbins(const TH1 *h, double epsilon)
Definition Resplot.cxx:838
simple struct to hold a value and it's associated uncertainty.
Definition StatVal.h:26
int r
Definition globals.cxx:22