KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
KVHistoManipulator.cpp
1//Created by KVClassFactory on Thu Oct 18 11:48:18 2007
2//Author: Eric Bonnet
3
4#include "KVHistoManipulator.h"
5#include "Riostream.h"
6#include "TProfile.h"
7#include "TList.h"
8#include "TString.h"
9#include "TClass.h"
10#include "TGraph.h"
11#include "TROOT.h"
12#include "TMethodCall.h"
13#include "TMath.h"
14#include "TRandom.h"
15#include "KVNumberList.h"
16#include "TSpline.h"
17#include "TStyle.h"
18#include "TCanvas.h"
19#include "TMultiGraph.h"
20
21#include <TObjString.h>
22
23using namespace std;
24
26
27KVHistoManipulator* gHistoManipulator;
28
29
32
34{
35 //Default constructor
36 init();
37 gHistoManipulator = this;
38 fVDCanvas = nullptr;
39}
40
41
42
43
45
47{
48 if (gHistoManipulator == this)
49 gHistoManipulator = nullptr;
50}
51
52
53//###############################################################################################################"
54//-------------------------------------------------
55
68
70{
71//-------------------------------------------------
72 // IMPORTANT l'histo est modifie
73 //
74 // Efface les bins dont la statistique est hors de l 'intervalle ]stat_min,stat_max[
75 // l'erreur et le contenu sont mis a zero
76 // le nombre d entree (GetEntries) reste inchange
77 // pour les TH1 cela correspond a GetBinContent(xx)
78 // pour les TH2 cela correspond a GetBinContent(xx,yy) --> cellules
79 // pour les TProfile cela correspond a GetBinEntries(xx)
80 // la fonction renvoie le nbre de bins ou cellules mis a zero
81 // si stat_min ou stat_max sont egales à -1 la borne associée n'est pas prise en compte
82 if (!hh) {
83 cout << "pointeur histogramme nul" << endl;
84 return -1;
85 }
86 Int_t nb_raz = 0;
87 Bool_t raz = kFALSE;
88 if (hh->InheritsFrom("TH2")) {
89
90 for (Int_t nx = 1; nx <= hh->GetNbinsX(); nx += 1) {
91 for (Int_t ny = 1; ny <= hh->GetNbinsY(); ny += 1) {
92 raz = kFALSE;
93 if (stat_min != -1 && hh->GetBinContent(nx, ny) < stat_min) raz = kTRUE;
94 if (stat_max != -1 && hh->GetBinContent(nx, ny) > stat_max) raz = kTRUE;
95 if (raz) {
96 hh->SetBinContent(nx, ny, 0);
97 hh->SetBinError(nx, ny, 0);
98 nb_raz += 1;
99 }
100 }
101 }
102 }
103 else if (hh->InheritsFrom("TProfile")) {
104 TProfile* prof = (TProfile*)hh;
105 for (Int_t nx = 1; nx <= prof->GetNbinsX(); nx += 1) {
106 raz = kFALSE;
107 if (stat_min != -1 && prof->GetBinEntries(nx) < stat_min) raz = kTRUE;
108 if (stat_max != -1 && prof->GetBinEntries(nx) > stat_max) raz = kTRUE;
109 if (raz) {
110 prof->SetBinContent(nx, 0);
111 prof->SetBinEntries(nx, 0);
112 prof->SetBinError(nx, 0);
113 nb_raz += 1;
114 }
115 }
116 }
117 else if (hh->InheritsFrom("TH1")) {
118 for (Int_t nx = 1; nx <= hh->GetNbinsX(); nx += 1) {
119 raz = kFALSE;
120 if (stat_min != -1 && hh->GetBinContent(nx) < stat_min) raz = kTRUE;
121 if (stat_max != -1 && hh->GetBinContent(nx) > stat_max) raz = kTRUE;
122 if (raz) {
123 hh->SetBinContent(nx, 0);
124 hh->SetBinError(nx, 0);
125 nb_raz += 1;
126 }
127 }
128 }
129
130 return nb_raz;
131}
132
133
134//###############################################################################################################"
135//-------------------------------------------------
136
146
148{
149//-------------------------------------------------
150 // IMPORTANT l'histo est modifie
151 //
152 // Applique une selection sur un histo 2D a partir d'un TCutG;
153 // Si mode=In seules les cellules comprises dans le TCutG sont gardes
154 // mode=Out --> Inverse
155 // la fonction renvoie le nbre de cellules mis a zero
156 // Attention le test ne se fait que sur les valeurs centrales des cellules (GetBinCenter())
157
158 if (!hh) {
159 cout << "pointeur histogramme nul" << endl;
160 return -1;
161 }
162 Int_t nb_raz = 0;
163 Bool_t raz = kFALSE;
164
165 for (Int_t nx = 1; nx <= hh->GetNbinsX(); nx += 1) {
166 for (Int_t ny = 1; ny <= hh->GetNbinsY(); ny += 1) {
167 raz = kFALSE;
168 Double_t valx = hh->GetXaxis()->GetBinCenter(nx);
169 Double_t valy = hh->GetYaxis()->GetBinCenter(ny);
170 if (mode == "in" && !cut->IsInside(valx, valy)) raz = kTRUE;
171 if (mode == "out" && cut->IsInside(valx, valy)) raz = kTRUE;
172 if (raz) {
173 hh->SetBinContent(nx, ny, 0);
174 hh->SetBinError(nx, ny, 0);
175 nb_raz += 1;
176 }
177 }
178 }
179
180 return nb_raz;
181}
182
183
184//###############################################################################################################
185//-------------------------------------------------
186
212
213TH1* KVHistoManipulator::ScaleHisto(TH1* hh, TF1* fx, TF1* fy, Int_t nx, Int_t ny, Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Option_t* norm)
214{
215//-------------------------------------------------
216 // Applique les transformations correspondantes aux fonctions (TF1 *) donnees en argument
217 // et donne en sortie l histogramme transforme celui ci a pour nom "nom_de_histo_input"_scaled
218 //
219 // IMPORTANT l'histo d'entree n'est pas modifie
220 //
221 // pour TH1 et TProfile n'est pris en compte que TF1 *fx;
222 // les parametres de la fonction doivent etre initialises avant.
223 // Si la fonction est un pointeur NULL, aucune transformation n est appliquee et l axe reste tel quel.
224 // L'intervalle de validite des fonctions TF1::SetRange est determine a partir des limites de l histo apres transformations
225 // nx xmin et xmax sont les parametres de binnage des histo
226 // si nx=-1, ceux ci sont recalcules automatiquement nx = nx_initial, xmin = fx(xmin), xmax = fx(xmax)
227 // si nx!=-1, xmin et xmax doivent etre aussi specifies
228 // il en est de meme pour l'axe y
229 //
230 // OPTION: norm
231 // norm = "" (default) : no adjustment is made for the change in bin width due to the transformation
232 // norm = "width" : bin contents are adjusted for width change, so that the integral of the histogram
233 // contents taking into account the bin width (i.e. TH1::Integral("width")) is the same.
234 //
235 // NOTE for 1D histos:
236 // If the binning of the scaled histogram is imposed by the user (nx!=-1), then in order
237 // to achieve a continuous scaled distribution we have to randomize X within each bin.
238 // This will only work if the bin contents of 'hh' are integer numbers, i.e. unnormalised raw data histogram.
239
240 Bool_t width = !strcmp(norm, "width");
241 TH1* gg = nullptr;
243 Bool_t fixed_bins = (nx != -1);
244 if (fx) {
245 if (nx == -1) {
246 nx = hh->GetNbinsX();
247 abs = hh->GetXaxis()->GetBinLowEdge(1);
248 xmin = fx->Eval(abs);
249 abs = hh->GetXaxis()->GetBinUpEdge(hh->GetNbinsX());
250 xmax = fx->Eval(abs);
251 }
252 fx->SetRange(hh->GetXaxis()->GetBinLowEdge(1), hh->GetXaxis()->GetBinUpEdge(nx));
253 fx->SetNpx(nx + 1);
254 }
255 else {
256 nx = hh->GetNbinsX();
257 xmin = hh->GetXaxis()->GetBinLowEdge(1);
258 xmax = hh->GetXaxis()->GetBinUpEdge(hh->GetNbinsX());
259 }
260
261 if (hh->InheritsFrom("TH2")) {
262 if (fy) {
263 if (ny == -1) {
264 ny = hh->GetNbinsY();
265 abs = hh->GetYaxis()->GetBinLowEdge(1);
266 ymin = fy->Eval(abs);
267 abs = hh->GetYaxis()->GetBinUpEdge(hh->GetNbinsY());
268 ymax = fy->Eval(abs);
269 }
270 fy->SetRange(hh->GetYaxis()->GetBinLowEdge(1), hh->GetYaxis()->GetBinUpEdge(ny));
271 fy->SetNpx(ny + 1);
272 }
273 else {
274 ny = hh->GetNbinsY();
275 ymin = hh->GetYaxis()->GetBinLowEdge(1);
276 ymax = hh->GetYaxis()->GetBinUpEdge(hh->GetNbinsY());
277 }
278 }
279
280 TClass* clas = TClass::GetClass(hh->ClassName());
281 gg = (TH1*) clas->New();
282 if (!gg) return nullptr;
283 TString hname;
284 hname.Form("%s_scaled", hh->GetName());
285 gg->SetNameTitle(hname.Data(), hh->GetTitle());
286
287 if (hh->InheritsFrom("TH2")) gg->SetBins(nx, xmin, xmax, ny, ymin, ymax);
288 else gg->SetBins(nx, xmin, xmax);
289
290 // if option norm="width', take account of change of X bin width between original and scaled histograms
291 Double_t Xbin_width_corr = 1.0, Ybin_width_corr = 1.0;
292 if (width) {
293 Double_t orig_Xbin_width = (hh->GetXaxis()->GetXmax() - hh->GetXaxis()->GetXmin()) / hh->GetNbinsX();
294 Double_t new_Xbin_width = (gg->GetXaxis()->GetXmax() - gg->GetXaxis()->GetXmin()) / gg->GetNbinsX();
295 Xbin_width_corr = orig_Xbin_width / new_Xbin_width;
296 if (hh->InheritsFrom("TH2")) {
297 // Take account of change of X bin width between original and scaled histograms
298 Double_t orig_Ybin_width = (hh->GetYaxis()->GetXmax() - hh->GetYaxis()->GetXmin()) / hh->GetNbinsY();
299 Double_t new_Ybin_width = (gg->GetYaxis()->GetXmax() - gg->GetYaxis()->GetXmin()) / gg->GetNbinsY();
300 Ybin_width_corr = orig_Ybin_width / new_Ybin_width;
301 }
302 }
303
304 for (Int_t xx = 1; xx <= hh->GetNbinsX(); xx += 1) {
305 Double_t bmin = hh->GetXaxis()->GetBinLowEdge(xx);
306 Double_t bmax = hh->GetXaxis()->GetBinUpEdge(xx);
307 abs = gRandom->Uniform(bmin, bmax);
308 if (abs == bmax) abs = bmin;
309 Double_t resx = abs;
310 if (fx) resx = fx->Eval(abs);
311 if (hh->InheritsFrom("TH2")) {
312 for (Int_t yy = 1; yy <= hh->GetNbinsY(); yy += 1) {
313 if (hh->GetBinContent(xx, yy) > 0) {
314 bmin = hh->GetYaxis()->GetBinLowEdge(yy);
315 bmax = hh->GetYaxis()->GetBinUpEdge(yy);
316 abs = gRandom->Uniform(bmin, bmax);
317 if (abs == bmax) abs = bmin;
318 Double_t resy = abs;
319 if (fy) resy = fy->Eval(abs);
320 gg->SetBinContent(gg->GetXaxis()->FindBin(resx),
321 gg->GetYaxis()->FindBin(resy),
322 hh->GetBinContent(xx, yy)*Xbin_width_corr * Ybin_width_corr);
323 //gg->SetBinError(gg->GetXaxis()->FindBin(resx),gg->GetYaxis()->FindBin(resy),hh->GetBinError(xx,yy));
324 }
325 }
326 }
327 else {
328 // 1-D histogram
329 if (fixed_bins) {
330 // histogram binning imposed by user. we need to randomise inside bins of original histogram
331 // otherwise scaled histogram will be discontinuously filled.
332 Int_t nmax = (Int_t)hh->GetBinContent(xx);
333 for (int i = 0; i < nmax; i++) {
334 abs = gRandom->Uniform(bmin, bmax);
335 Double_t resx = fx->Eval(abs);
336 gg->Fill(resx, Xbin_width_corr);
337 //cout << resx << " " << Xbin_width_corr << endl;
338 }
339 }
340 else {
341 // range and binning calculated from function.
342 Double_t resy = hh->GetBinContent(xx);
343 if (fy) resy = fy->Eval(resy);
344 gg->SetBinContent(gg->GetXaxis()->FindBin(resx), resy * Xbin_width_corr);
345 }
346 }
347 }
348 return gg;
349}
350
351
352
359
360TGraph* KVHistoManipulator::ScaleGraph(const TGraph* hh, const TString& axis, const TF1& f1, const TF1& f2) const
361{
362 // Returns pointer to new graph whose X- and/or Y-coordinates are those of the input graph scaled according to the given function(s).
363 //
364 // * If axis="X", f1 is applied to X-coordinates
365 // * If axis="Y", f1 is applied to Y-coordinates
366 // * If axis="XY", f1 is applied to X-coordinates and f2 is applied to Y-coordinates
367
368 TGraph* gg = nullptr;
369 TClass* clas = TClass::GetClass(hh->ClassName());
370 gg = (TGraph*) clas->New();
371 if (!gg) return nullptr;
372 TString hname;
373 hname.Form("%s_scaled", hh->GetName());
374 gg->SetNameTitle(hname.Data(), hh->GetTitle());
375
376 bool scale_X = axis.Contains("X");
377 const TF1* FX = &f1;
378 bool scale_Y = axis.Contains("Y");
379 const TF1* FY = &f2;
380 if (axis == "Y") FY = &f1;
381
382 Int_t np = hh->GetN();
383 for (Int_t nn = 0; nn < np; nn += 1) {
384 Double_t xx1, yy1;
385 hh->GetPoint(nn, xx1, yy1);
386 Double_t xx2 = xx1;
387 if (scale_X) xx2 = FX->Eval(xx1);
388 Double_t yy2 = yy1;
389 if (scale_Y) yy2 = FY->Eval(yy1);
390 gg->SetPoint(nn, xx2, yy2);
391 if (gg->InheritsFrom("TGraphErrors")) {
392 // transformation of errors
393 // if f = f(x) the error on f, e_f, for a given error on x, e_x, is
394 // e_f = abs0(df/dx) * e_x
395 Double_t e_x = ((TGraphErrors*)hh)->GetErrorX(nn);
396 Double_t e_y = ((TGraphErrors*)hh)->GetErrorY(nn);
397 if (scale_X) e_x = TMath::Abs(FX->Derivative(xx1)) * e_x;
398 if (scale_Y) e_y = TMath::Abs(FY->Derivative(yy1)) * e_y;
399 ((TGraphErrors*)gg)->SetPointError(nn, e_x, e_y);
400 }
401 }
402
403 return gg;
404}
405
406//###############################################################################################################
407//-------------------------------------------------
408
413
415{
416//-------------------------------------------------
417 //Exemple d utilisation de la methode KVHistoManipulator::ScaleHisto avec ici
418 //L'obtention des distributions centrees reduites
419
420 if (!hh) {
421 cout << "pointeur histogramme nul" << endl;
422 return hh;
423 }
424
425 TString expression;
426 expression.Form("(x-%lf)/%lf", hh->GetMean(1), hh->GetRMS(1));
427 TF1 fx("fx", expression.Data());
428 unique_ptr<TF1> fy;
429 if (hh->InheritsFrom("TH2")) {
430 expression.Form("(x-%lf)/%lf", hh->GetMean(2), hh->GetRMS(2));
431 fy.reset(new TF1("fy", expression.Data()));
432 }
433
434 TH1* gg = ScaleHisto(hh, &fx, fy.get(), nx, ny, xmin, xmax, ymin, ymax);
435 TString hname;
436 hname.Form("%s_centred", hh->GetName());
437 gg->SetName(hname.Data());
438
439 return gg;
440}
441
442//###############################################################################################################
443//-------------------------------------------------
444
449
451{
452//-------------------------------------------------
453 //Exemple d utilisation de la methode KVHistoManipulator::ScaleHisto avec ici
454 //L'obtention des distributions centrees reduites
455
456 if (!hh) {
457 cout << "pointeur histogramme nul" << endl;
458 return hh;
459 }
460
461 TString expression;
462 expression.Form("(x-%lf)/%lf", hh->GetMean(1), hh->GetRMS(1));
463 TF1 fx("fx", expression.Data());
464 TH2* gg = (TH2*)ScaleHisto(hh, &fx, nullptr, nx, -1, xmin, xmax, -1., -1.);
465 TString hname;
466 hname.Form("%s_centred_X", hh->GetName());
467 gg->SetName(hname.Data());
468
469 return gg;
470}
471
472//###############################################################################################################
473//-------------------------------------------------
474
479
481{
482//-------------------------------------------------
483 //Exemple d utilisation de la methode KVHistoManipulator::ScaleHisto avec ici
484 //L'obtention des distributions centrees reduites
485
486 if (!hh) {
487 cout << "pointeur histogramme nul" << endl;
488 return hh;
489 }
490
491 TString expression;
492 expression.Form("(x-%lf)/%lf", hh->GetMean(1), hh->GetRMS(1));
493 TF1 fy("fy", expression.Data());
494 TH2* gg = (TH2*)ScaleHisto(hh, nullptr, &fy, -1, ny, -1., -1., ymin, ymax);
495 TString hname;
496 hname.Form("%s_centred_Y", hh->GetName());
497 gg->SetName(hname.Data());
498
499 return gg;
500}
501
502
503//###############################################################################################################
504//-------------------------------------------------
505
520
522{
523//-------------------------------------------------
524 // Renormalisation de l histogramme 2D sous la contrainte d'une ditribution plate suivant X ou Y
525 // (signal de bimodalite par exemple)
526 // et donne en sortie l histogramme transforme celui ci a pour nom "nom_de_histo_input"_normalised
527 //
528 // IMPORTANT l'histo d'entree n'est pas modifie
529 //
530 // bmin, bmax : intervalle de bins ou l'on effectue la renormalisation (par defaut -1,-1 --> toute la plage)
531 // les contenus hors de l intervalle [bmin,bmax] sont remis a zero
532 // Si on veut une distribution plate en X ils indiquent l'intervalle sur cet axe x
533 // valref : valeur en stat de la distribution plate (par defaut 1)
534 //
535 // ATTENTION la propagation des erreurs n est pas (pour l instant) implementee
536
537 if (!hh) {
538 cout << "pointeur histogramme nul" << endl;
539 return hh;
540 }
541
542 if (bmin == -1) bmin = 1;
543 TString hname;
544 hname.Form("%s_normalised", hh->GetName());
545 TH2* clone = 0;
546 if ((clone = (TH2*)gDirectory->Get(hname.Data()))) delete clone;
547 clone = (TH2*)hh->Clone(hname.Data());
548 clone->Reset();
549
550 if (axis == "X") {
551 if (bmax == -1) {
552 bmax = hh->GetNbinsX();
553 }
554 for (Int_t nx = bmin; nx <= bmax; nx += 1) {
555 Double_t integ = 0;
556 for (Int_t ny = 1; ny <= hh->GetNbinsY(); ny += 1) integ += hh->GetBinContent(nx, ny);
557 if (integ > 0) {
558 for (Int_t ny = 1; ny <= hh->GetNbinsY(); ny += 1) {
559 clone->SetBinContent(nx, ny, hh->GetBinContent(nx, ny)*valref / integ);
560 if (hh->GetBinContent(nx, ny) > 0) {
561 Double_t erreur = clone->GetBinContent(nx, ny) * TMath::Sqrt(1. / hh->GetBinContent(nx, ny) + 1. / integ);
562 clone->SetBinError(nx, ny, erreur);
563 }
564 }
565 }
566 }
567 }
568 else if (axis == "Y") {
569 if (bmax == -1) {
570 bmax = hh->GetNbinsY();
571 }
572 for (Int_t ny = bmin; ny <= bmax; ny += 1) {
573 Double_t integ = 0;
574 for (Int_t nx = 1; nx <= hh->GetNbinsX(); nx += 1) integ += hh->GetBinContent(nx, ny);
575 if (integ > 0) {
576 for (Int_t nx = 1; nx <= hh->GetNbinsX(); nx += 1) {
577 clone->SetBinContent(nx, ny, hh->GetBinContent(nx, ny)*valref / integ);
578 Double_t erreur = clone->GetBinContent(nx, ny) * TMath::Sqrt(1. / hh->GetBinContent(nx, ny) + 1. / integ);
579 clone->SetBinError(nx, ny, erreur);
580 }
581 }
582 }
583 }
584 else {
585 cout << "l option TString axis doit etre X ou Y" << endl;
586 }
587 return clone;
588}
589
590//###############################################################################################################
591//-------------------------------------------------
592
598
600{
601//-------------------------------------------------
602 // Les bornes valmin, valmax definissent l'intervalle ou l on effectue la renormalisation
603 // On en derive les valeurs de bins
604 // pour la methode KVHistoManipulator::RenormaliseHisto(TH1 *hh,TString axis,Double_t valref,Int_t bmin,Int_t bmax)
605
606 if (!hh) {
607 cout << "pointeur histogramme nul" << endl;
608 return hh;
609 }
610 if (axis == "X") return RenormaliseHisto(hh, hh->GetXaxis()->FindBin(valmin), hh->GetXaxis()->FindBin(valmax), axis, valref);
611 else if (axis == "Y") return RenormaliseHisto(hh, hh->GetYaxis()->FindBin(valmin), hh->GetYaxis()->FindBin(valmax), axis, valref);
612 else {
613 cout << "l option TString axis doit etre X ou Y" << endl;
614 return nullptr;
615 }
616
617}
618
619
620//###############################################################################################################
621//-------------------------------------------------
622
640
642{
643//-------------------------------------------------
644 // Cumule le contenu de l histo hh entre bin bmin et bmax et retourne l histo correspondant
645 // Si direction="C" (default):
646 // SetBinContent(n) = GetBinContent(bmin)+GetBinContent(bmin+1)+ ... + GetBinContent(n)
647 // Si direction="D":
648 // SetBinContent(n) = GetBinContent(bmax)+GetBinContent(bmax-1)+ ... + GetBinContent(n)
649 //
650 // Donne en sortie l histogramme transforme celui ci a pour nom "nom_de_histo_input"_cumulated
651 //
652 // IMPORTANT l'histo d'entree n'est pas modifie
653 //
654 // si bmin=-1, bmin=1
655 // si bmax=-1, bmax=GetNbinsX()
656 //
657 // Avec norm = "surf" (default) l'integral de l histo cumule est egale a 1
658 // Avec norm = "max" le contenu de l'histogram est renormalise de facon a ce que le maximum soit 1
659
660 if (!hh) {
661 cout << "pointeur histogramme nul" << endl;
662 return hh;
663 }
664 direction.ToUpper();
665 if (direction != "C" && direction != "D") {
666 cout << "l option TString direction doit etre C ou D" << endl;
667 return nullptr;
668 }
669 if (hh->GetDimension() == 1) {
670 if (bmin < 1) bmin = 1;
671 if (bmax < 1) bmax = hh->GetNbinsX();
672 TString hname;
673 hname.Form("%s_cumulated", hh->GetName());
674 TH1* clone = (TH1*)hh->Clone(hname.Data());
675 clone->Reset();
676 Double_t sum = 0;
677 Double_t big_sum = 0;
678 if (direction == "C") {
679 /* "standard" cumulative histogram, from bin 'bmin' to bin 'bmax' */
680 for (Int_t nx = 1; nx <= hh->GetNbinsX(); ++nx) {
681 if (nx < bmin) clone->SetBinContent(nx, 0);
682 else if (nx > bmax) clone->SetBinContent(nx, sum);
683 else {
684 sum += hh->GetBinContent(nx);
685 clone->SetBinContent(nx, sum);
686 }
687 big_sum += sum;
688 }
689 }
690 else {
691 /* "reverse" cumulative histogram, from bin 'bmax' to bin 'bmin' */
692 for (Int_t nx = hh->GetNbinsX(); nx >= 1; --nx) {
693 if (nx > bmax) clone->SetBinContent(nx, 0);
694 else if (nx < bmin) clone->SetBinContent(nx, sum);
695 else {
696 sum += hh->GetBinContent(nx);
697 clone->SetBinContent(nx, sum);
698 }
699 big_sum += sum;
700 }
701 }
702 // normalisation
703 if (!strcmp(norm, "surf")) {
704 clone->Scale(1. / big_sum);
705 }
706 else if (!strcmp(norm, "max")) {
707 clone->Scale(1. / sum);
708 }
709 return clone;
710 }
711 else {
712 cout << "cette methode n est pas prevue pour les TH2, TH3" << endl;
713 return nullptr;
714 }
715
716}
717
718
719//###############################################################################################################
720//-------------------------------------------------
721
734
736{
737//-------------------------------------------------
738 // retourne la derivee d ordre 0, 1 ou 2 d'un histogramme
739 // 0 -> correspond a un lissage (smooth) de la distribution
740 // 1 et 2 correspondent aux derivees premieres et deuxiemes
741 //
742 // 0 -> derivee zero Yi = 1/35*(-3yi-2 + 12yi-1 +17yi +12yi1 -3yi2)
743 // 1 -> derivee premiere Y'i = 1/12h*(yi-2 - 8yi-1 +8yi1 -1yi2)
744 // 2 -> derivee seconde Y''i = 1/7h/h*(2yi-2 -1 yi-1 -2yi -1yi1 +2yi2)
745 // les derivees pour les bins 1,2 et n-1,n ne sont pas calculees
746 //
747 // IMPORTANT l'histo d'entree n'est pas modifie
748
749 if (!hh) {
750 cout << "pointeur histogramme nul" << endl;
751 return hh;
752 }
753 if (!(0 <= order && order <= 2)) {
754 cout << "ordre " << order << "n est pas implemente" << endl;
755 return nullptr;
756 }
757 if (hh->GetDimension() == 1) {
758
759 TString hname;
760 hname.Form("%s_derivated_%d", hh->GetName(), order);
761 TH1* clone = 0;
762
763 if (hh->InheritsFrom("TProfile")) clone = new TH1F(hname.Data(), hh->GetTitle(), hh->GetNbinsX(), hh->GetBinLowEdge(1), hh->GetBinLowEdge(hh->GetNbinsX() + 1));
764 else clone = (TH1*)hh->Clone(hname.Data());
765
766 clone->Reset();
767
768 for (Int_t nx = 3; nx <= hh->GetNbinsX() - 2; nx += 1) {
769 Double_t dev = 0;
770 if (order == 0) {
771 dev = 1. / 35.*(
772 -3 * hh->GetBinContent(nx - 2)
773 + 12 * hh->GetBinContent(nx - 1)
774 + 17 * hh->GetBinContent(nx)
775 + 12 * hh->GetBinContent(nx + 1)
776 - 3 * hh->GetBinContent(nx + 2)
777 );
778 }
779 else if (order == 1) {
780 Double_t h = hh->GetBinWidth(1);
781 dev = 1 / 12. / h * (
782 1 * hh->GetBinContent(nx - 2)
783 - 8 * hh->GetBinContent(nx - 1)
784 + 0 * hh->GetBinContent(nx)
785 + 8 * hh->GetBinContent(nx + 1)
786 - 1 * hh->GetBinContent(nx + 2)
787 );
788 }
789 else {
790 Double_t h2 = pow(hh->GetBinWidth(1), 2.);
791 dev = 1 / 7. / h2 * (
792 2 * hh->GetBinContent(nx - 2)
793 - 1 * hh->GetBinContent(nx - 1)
794 - 2 * hh->GetBinContent(nx)
795 - 1 * hh->GetBinContent(nx + 1)
796 + 2 * hh->GetBinContent(nx + 2)
797 );
798 }
799 clone->SetBinContent(nx, dev);
800 }
801 return clone;
802 }
803 else {
804 cout << "cette methode n est pas prevue pour les TH2, TH3" << endl;
805 return nullptr;
806 }
807
808}
809
810
811//###############################################################################################################
812//-------------------------------------------------
813
830
832{
833//-------------------------------------------------
834 // Renvoie un TGraph
835 // A partir d'un 2D histo, permet de tracer l'evolution d'un moment en fonction d'un autre moment
836 // ou d'un moment d'un axe en fonction de l'aure axe.
837 // Les TString momentx et momenty doivent etre des "Get like" methodes connues de TH1
838 // comme GetMean, GetRMS ou GetKurtosis :
839 // - Si momenty="" -> on obtient l'evolution du moment momentx en fonction de
840 // la variable associée a l'autre axe.
841 // - Si momenty!="" -> on obtient l'evolution du moment momenty en fonction du momentx
842 //
843 // Si TString axis = X la variable en X est le parametre d echantillonage et vice-versa
844 //
845 // Exemple :
846 // GetMomentEvolution(histo,"GetMean","GetSkewness","X") -> Evolution de la skewness en fonction de la moyenne de l'observable
847 // en Y par tranche de l'observable X
848
849 if (axis != "X" && axis != "Y") {
850 cout << "GetMomentEvolution(TH2*,TString ,TString ,TString) Mauvaise syntaxe pour TString axis (X ou Y) " << endl;
851 return nullptr;
852 }
853 TMethodCall cmx;
854 cmx.InitWithPrototype(TClass::GetClass("TH1D"), Form("%s", momentx.Data()), "int");
855 if (!cmx.IsValid()) {
856 cout << "GetMomentEvolution(TH2*,TString ,TString ,TString) TString momentx n'est pas une methode valide " << momentx.Data() << endl;
857 return nullptr;
858 }
859 unique_ptr<TMethodCall> Ecmx(new TMethodCall());
860 Ecmx->InitWithPrototype(TClass::GetClass("TH1D"), Form("%sError", momentx.Data()), "int");
861
862 unique_ptr<TMethodCall> cmy, Ecmy;
863 if (momenty != "") {
864 cmy.reset(new TMethodCall());
865 cmy->InitWithPrototype(TClass::GetClass("TH1D"), Form("%s", momenty.Data()), "int");
866 if (!cmy->IsValid()) {
867 cout << "GetMomentEvolution(TH2*,TString ,TString ,TString) TString momenty n'est pas une methode valide " << momenty.Data() << endl;
868 return nullptr;
869 }
870 Ecmy.reset(new TMethodCall());
871 Ecmy->InitWithPrototype(TClass::GetClass("TH1D"), Form("%sError", momenty.Data()), "int");
872 }
873
874 TString fmt_histo;
875 fmt_histo.Form("GetMomentEvolution_%s", hh->GetName());
876 KVNumberList lbins = "";
877 Int_t nmax = 0;
878 if (axis == "Y") nmax = hh->GetNbinsX();
879 else nmax = hh->GetNbinsY();
880
881 Int_t npts = 0;
882 for (Int_t nn = 1; nn <= nmax; nn += 1) {
883 Double_t stat = 0;
884 if (axis == "Y") stat = ((TH1D*)hh->ProjectionY(fmt_histo.Data(), nn, nn))->Integral();
885 else stat = ((TH1D*)hh->ProjectionX(fmt_histo.Data(), nn, nn))->Integral();
886 if (stat > stat_min) {
887 npts += 1;
888 lbins.Add(nn);
889 }
890 gDirectory->Delete(fmt_histo.Data());
891 }
892 if (npts > 0) {
893 TGraphErrors* gr = new TGraphErrors(npts);
894 TH1D* hp;
895 npts = 0;
896 Double_t valx, valy, Evaly = 0, Evalx = 0;
897 lbins.Begin();
898 while (!lbins.End()) {
899 Int_t bin = lbins.Next();
900 if (axis == "Y") hp = (TH1D*)hh->ProjectionY(fmt_histo.Data(), bin, bin);
901 else hp = (TH1D*)hh->ProjectionX(fmt_histo.Data(), bin, bin);
902 cmx.Execute(hp, "1", valx);
903 if (Ecmx->IsValid()) Ecmx->Execute(hp, "1", Evalx);
904 if (momenty != "") {
905 cmy->Execute(hp, "1", valy);
906 if (Ecmy->IsValid()) Ecmy->Execute(hp, "1", Evaly);
907 gr->SetPoint(npts, valx, valy);
908 if (Evalx != 0 && Evaly != 0) gr->SetPointError(npts, Evalx, Evaly);
909 npts += 1;
910 }
911 else {
912 if (axis == "Y") {
913 valy = hh->GetXaxis()->GetBinCenter(bin);
914 Evaly = hh->GetXaxis()->GetBinWidth(bin) / 2.;
915 }
916 else {
917 valy = hh->GetYaxis()->GetBinCenter(bin);
918 Evaly = hh->GetYaxis()->GetBinWidth(bin) / 2.;
919 }
920 gr->SetPoint(npts, valy, valx);
921 if (Evalx != 0) gr->SetPointError(npts, Evaly, Evalx);
922 npts += 1;
923 }
924 gDirectory->Delete(fmt_histo.Data());
925 }
926 return gr;
927 }
928 else {
929 cout << "GetMomentEvolution(TH2*,TString ,TString ,TString) Aucun point dans le TGraph" << endl;
930 return nullptr;
931 }
932
933}
934
935
936
937//-------------------------------------------------
938
945
947{
948//-------------------------------------------------
949 // A partir de deux graphs ayant le meme nombre de points et le meme axe X,
950 // cette methode produit un graph correspondant a la correlation entre les deux axes Y
951 // Les inputs peuvent etre aussi bien des TGraph que des TGraphErrors dans ce dernier cas
952 // les barres d erreurs sont prises en compte
953
954 Double_t* yy = gry->GetY();
955 Double_t* xx = grx->GetY();
956 Double_t* ey = 0;
957 Double_t* ex = 0;
958
959 if (grx->GetN() != gry->GetN()) {
960 printf("ERREUR : KVHistoManipulator::LinkGraphs : les deux graphs n ont pas le meme nbre de points\n");
961 return 0;
962 }
963 Int_t npoints = grx->GetN();
964
965 TGraph* corre = 0;
966 if (grx->InheritsFrom("TGraphErrors") || gry->InheritsFrom("TGraphErrors")) {
967 if (grx->InheritsFrom("TGraphErrors")) ex = grx->GetEY();
968 if (gry->InheritsFrom("TGraphErrors")) ey = gry->GetEY();
969 corre = new TGraphErrors(npoints, xx, yy, ex, ey);
970 }
971 else corre = new TGraph(npoints, xx, yy);
972
974 name.Form("corre_%s_VS_%s", gry->GetName(), grx->GetName());
975 corre->SetNameTitle(name.Data(), grx->GetTitle());
976
977 return corre;
978
979}
980
981
982
983//###############################################################################################################
984//-------------------------------------------------
985
990
992{
993//-------------------------------------------------
994 // Calcul le coefficient de correlation entre les variables X et Y
995 // d'un bidim... Equivalent a TH2F::GetCorrelationFactor() de ROOT
996 Double_t sumxy = 0;
997 Double_t sumx = 0, sumx2 = 0;
998 Double_t sumy = 0, sumy2 = 0;
999
1000 Double_t compt = 0;
1001 for (Int_t nx = 1; nx <= hh->GetNbinsX(); nx += 1) {
1002 for (Int_t ny = 1; ny <= hh->GetNbinsY(); ny += 1) {
1003 compt += hh->GetBinContent(nx, ny);
1004 sumxy += hh->GetBinContent(nx, ny) * hh->GetXaxis()->GetBinCenter(nx) * hh->GetYaxis()->GetBinCenter(ny);
1005 sumx += hh->GetBinContent(nx, ny) * hh->GetXaxis()->GetBinCenter(nx);
1006 sumy += hh->GetBinContent(nx, ny) * hh->GetYaxis()->GetBinCenter(ny);
1007 sumx2 += hh->GetBinContent(nx, ny) * pow(hh->GetXaxis()->GetBinCenter(nx), 2.);
1008 sumy2 += hh->GetBinContent(nx, ny) * pow(hh->GetYaxis()->GetBinCenter(ny), 2.);
1009 }
1010 }
1011
1012 Double_t meanxy = sumxy / compt;
1013 Double_t meanx = sumx / compt;
1014 Double_t meany = sumy / compt;
1015 Double_t meanx2 = sumx2 / compt;
1016 Double_t sigmax = sqrt(meanx2 - pow(meanx, 2.));
1017 Double_t meany2 = sumy2 / compt;
1018 Double_t sigmay = sqrt(meany2 - pow(meany, 2.));
1019
1020 Double_t rho = (meanxy - meanx * meany) / (sigmax * sigmay);
1021 return rho;
1022
1023}
1024
1025//###############################################################################################################"
1026//-------------------------------------------------
1027
1033
1035{
1036//-------------------------------------------------
1037 // Retourne une liste contenant les projections par tranche de l'axe (TString axis="X" ou "Y")
1038 // remplissant une condition leur integral qui doit etre superieur à MinIntegral (=-1 par defaut)
1039 // si axis="X", les projections sur l'axe Y de l'histogramme est fait pour chaque bin de l'axe X
1040
1041 if (!hh) {
1042 cout << "pointeur histogramme nul" << endl;
1043 return NULL;
1044 }
1045 TH1D* h1d = 0;
1046 TString proj_name;
1047 if (hh->InheritsFrom("TH2")) {
1048 KVList* list = new KVList();
1049 cout << "TH2" << endl;
1050 if (axis == "X") {
1051 for (Int_t nx = 1; nx <= hh->GetNbinsX(); nx += 1) {
1052 Double_t integ = 0;
1053 for (Int_t ny = 1; ny <= hh->GetNbinsY(); ny += 1)
1054 integ += hh->GetBinContent(nx, ny);
1055
1056 if (integ > MinIntegral) {
1057 proj_name.Form("%s_bX_%d", hh->GetName(), nx);
1058 h1d = hh->ProjectionY(proj_name.Data(), nx, nx);
1059 h1d->SetTitle(Form("%lf", hh->GetXaxis()->GetBinCenter(nx)));
1060 list->Add(h1d);
1061 }
1062
1063 }
1064 }
1065 else if (axis == "Y") {
1066 for (Int_t ny = 1; ny <= hh->GetNbinsY(); ny += 1) {
1067 Double_t integ = 0;
1068 for (Int_t nx = 1; nx <= hh->GetNbinsX(); nx += 1)
1069 integ += hh->GetBinContent(nx, ny);
1070
1071 if (integ > MinIntegral) {
1072 proj_name.Form("%s_bY_%d", hh->GetName(), ny);
1073 h1d = hh->ProjectionX(proj_name.Data(), ny, ny);
1074 h1d->SetTitle(Form("%lf", hh->GetYaxis()->GetBinCenter(ny)));
1075 list->Add(h1d);
1076 }
1077 }
1078 }
1079 else {
1080 cout << "l option TString axis doit etre X ou Y" << endl;
1081 }
1082
1083 return list;
1084 }
1085 else {
1086 cout << "cette methode n est prevue que pour les TH2 and sons" << endl;
1087 return NULL;
1088 }
1089
1090}
1091
1092
1093//-------------------------------------------------
1094
1103
1105{
1106//-------------------------------------------------
1107 // Donne les valeurs permettant de decouper l'integral
1108 // d'un histograme en tranche de 10% de la stat total
1109 // retourne une KVNumberList indiquant les delimitation des tranches
1110 //
1111 // L'utilisateur doit effacer cette liste apres utilisation
1112 //
1113
1114 if (!hh) {
1115 cout << "pointeur histogramme nul" << endl;
1116 return NULL;
1117 }
1118
1119 TString proj_name;
1120 Double_t integral = 0;
1121 for (Int_t nx = 1; nx <= hh->GetXaxis()->GetNbins(); nx += 1) {
1122 integral += hh->GetBinContent(nx);
1123 }
1124
1125 KVNumberList* nl = new KVNumberList();
1126
1127 Double_t tranche = 0;
1128 printf("integral du spectre %lf -> tranche de %lf\n", integral, integral / ntranches);
1129
1130 Int_t nt = 0;
1131 for (Int_t nx = hh->GetXaxis()->GetNbins(); nx >= 1; nx -= 1) {
1132 tranche += hh->GetBinContent(nx);
1133 //printf("nx=%d, nt=%d, tranche[nt]=%lf -> %lf\n",nx,nt,tranche,integral/ntranches);
1134 if (tranche >= integral / ntranches) {
1135
1136 nt += 1;
1137 //printf("\tnouvelle tranche %d %d\n",nx,nt);
1138 nl->Add(nx);
1139 tranche = 0;
1140
1141 }
1142 }
1143
1144 //Verif
1145 /*
1146 Int_t sup=1;
1147 Int_t inf = 1;
1148 nl->Begin();
1149 while (!nl->End()){
1150 sup = nl->Next();
1151 printf("%d %d - %lf -> %lf\n",inf,sup,hh->Integral(inf,sup),hh->Integral(inf,sup)/integral*100);
1152 inf = sup+1;
1153 }
1154 sup = hh->GetXaxis()->GetNbins();
1155 printf("%d %d - %lf -> %lf\n",inf,sup,hh->Integral(inf,sup),hh->Integral(inf,sup)/integral*100);
1156 */
1157 return nl;
1158}
1159
1160
1161//-------------------------------------------------
1162
1169
1171{
1172//-------------------------------------------------
1173 // Cree un histo 2D en intervertissant les axes
1174 //
1175 // L'utilisateur doit effacer cet histo apres utilisation
1176 //
1177
1178 if (!hh) {
1179 cout << "pointeur histogramme nul" << endl;
1180 return NULL;
1181 }
1182 if (!hh->InheritsFrom("TH2")) {
1183 Error("PermuteAxis", "methode definie uniquement pour les classes TH2 et filles");
1184 }
1185 Int_t nx = hh->GetNbinsX();
1186 Int_t ny = hh->GetNbinsY();
1187
1188 TH2F* gg = new TH2F(
1189 Form("%s_perm", hh->GetName()),
1190 hh->GetTitle(),
1191 ny,
1192 hh->GetYaxis()->GetBinLowEdge(1),
1193 hh->GetYaxis()->GetBinLowEdge(ny + 1),
1194 nx,
1195 hh->GetXaxis()->GetBinLowEdge(1),
1196 hh->GetXaxis()->GetBinLowEdge(nx + 1)
1197 );
1198
1199 for (Int_t xx = 1; xx <= nx; xx += 1) {
1200 for (Int_t yy = 1; yy <= ny; yy += 1) {
1201
1202 gg->SetBinContent(yy, xx, hh->GetBinContent(xx, yy));
1203
1204 }
1205 }
1206 return gg;
1207}
1208
1209
1210//-------------------------------------------------
1211
1218
1220{
1221//-------------------------------------------------
1222 // Cree un TGraph en intervertissant les axes
1223 //
1224 // L'utilisateur doit effacer ce graph apres utilisation
1225 //
1226
1227 if (!gr) {
1228 cout << "pointeur graph nul" << endl;
1229 return NULL;
1230 }
1231 if (!gr->InheritsFrom("TGraph")) {
1232 Error("PermuteAxis", "methode definie uniquement pour les classes TGraph et filles");
1233 }
1234
1235 TGraphErrors* gr2 = new TGraphErrors();
1236 for (Int_t nn = 0; nn < gr->GetN(); nn += 1) {
1237 Double_t px, py;
1238 gr->GetPoint(nn, px, py);
1239 gr2->SetPoint(nn, py, px);
1240 if (gr->InheritsFrom("TGraphErrors")) {
1241 gr2->SetPointError(nn, ((TGraphErrors*)gr)->GetErrorY(nn), ((TGraphErrors*)gr)->GetErrorX(nn));
1242 }
1243 }
1244
1245 return gr2;
1246}
1247
1248
1249//-------------------------------------------------
1250
1257
1259{
1260//-------------------------------------------------
1261 // Cree un graph à partir d un histo
1262 //
1263 // L'utilisateur doit effacer ce TGraph apres utilisation
1264 //
1265
1266 if (!pf) {
1267 cout << "pointeur histogramme nul" << endl;
1268 return NULL;
1269 }
1270 if (!pf->InheritsFrom("TProfile")) {
1271 //Error("MakeGraphFrom","methode definie uniquement pour les classes TProfile et filles");
1272 return NULL;
1273 }
1274 Int_t nx = pf->GetNbinsX();
1275
1276 TGraphErrors* gr = new TGraphErrors();
1277 for (Int_t xx = 1; xx <= nx; xx += 1) {
1278 if (pf->GetBinEntries(xx) > 0) {
1279 gr->SetPoint(gr->GetN(), pf->GetBinCenter(xx), pf->GetBinContent(xx));
1280 if (Error)
1281 gr->SetPointError(gr->GetN() - 1, pf->GetBinWidth(xx) / 2, pf->GetBinError(xx));
1282 }
1283 }
1284
1285 return gr;
1286}
1287
1288
1289
1301
1303{
1304 // Extract mean and standard deviation (if generated with mode "profs") or error on mean (with "prof")
1305 // from TProfile and put their evolution in 2 TGraphs.
1306 // \returns Pointer to graph containing mean values
1307 // \param sigma Pointer to graph containing standard deviation
1308 //
1309 // Example of use:
1310 //~~~~{.cpp}
1311 // TGraph* sig_graf;
1312 // TGraph* mean_graf = HM.ExtractMeanAndSigmaFromProfile(pf, sig_graf);
1313 //~~~~
1314
1315 Int_t nx = pf->GetNbinsX();
1316
1317 TGraph* gr = new TGraph;
1318 sigma = new TGraph;
1319 int i = 0;
1320 for (Int_t xx = 1; xx <= nx; xx += 1) {
1321 if (pf->GetBinEntries(xx) > 0) {
1322 gr->SetPoint(i, pf->GetBinCenter(xx), pf->GetBinContent(xx));
1323 sigma->SetPoint(i, pf->GetBinCenter(xx), pf->GetBinError(xx));
1324 ++i;
1325 }
1326 }
1327
1328 return gr;
1329}
1330
1331
1332// //-------------------------------------------------
1333// TGraphErrors* KVHistoManipulator::MakeGraphFrom(TH1* pf, Bool_t Error)
1334// {
1335// //-------------------------------------------------
1336// // Cree un graph à partir d un histo
1337// //
1338// // L'utilisateur doit effacer ce TGraph apres utilisation
1339// //
1340//
1341// if (!pf) {
1342// cout << "pointeur histogramme nul" << endl;
1343// return NULL;
1344// }
1345// if (!pf->InheritsFrom("TH1")) {
1346// //Error("MakeGraphFrom","methode definie uniquement pour les classes TProfile et filles");
1347// }
1348// else if (pf->InheritsFrom("TProfile"))
1349// {
1350// return
1351// }
1352// Int_t nx = pf->GetNbinsX();
1353//
1354// TGraphErrors* gr = new TGraphErrors();
1355// for (Int_t xx = 1; xx <= nx; xx += 1) {
1356// if (pf->GetBinEntries(xx) > 0) {
1357// gr->SetPoint(gr->GetN(), pf->GetBinCenter(xx), pf->GetBinContent(xx));
1358// if (Error)
1359// gr->SetPointError(gr->GetN() - 1, pf->GetBinWidth(xx) / 2, pf->GetBinError(xx));
1360// }
1361// }
1362//
1363// return gr;
1364// }
1365
1366//###############################################################################################################"
1367//-------------------------------------------------
1368
1371
1372void KVHistoManipulator::DefinePattern(TH1* ob, TString titleX, TString titleY, TString labelX, TString labelY)
1373{
1374//-------------------------------------------------
1375 DefinePattern(ob->GetXaxis(), titleX, labelX);
1376 DefinePattern(ob->GetYaxis(), titleY, labelY);
1377
1378}
1379
1380//###############################################################################################################"
1381//-------------------------------------------------
1382
1385
1387{
1388//-------------------------------------------------
1389 DefinePattern(ob->GetXaxis(), titleX, labelX);
1390 DefinePattern(ob->GetYaxis(), titleY, labelY);
1391
1392}
1393
1394//###############################################################################################################"
1395//-------------------------------------------------
1396
1399
1400void KVHistoManipulator::DefinePattern(TF1* ob, TString titleX, TString titleY, TString labelX, TString labelY)
1401{
1402//-------------------------------------------------
1403 DefinePattern(ob->GetXaxis(), titleX, labelX);
1404 DefinePattern(ob->GetYaxis(), titleY, labelY);
1405
1406}
1407
1408//###############################################################################################################"
1409//-------------------------------------------------
1410
1413
1415{
1416//-------------------------------------------------
1417
1418 TObjArray* tok = NULL;
1419
1420 if (!title.IsNull()) {
1421 tok = title.Tokenize(" ");
1422 if (tok->GetEntries() == 3) {
1423 ax->SetTitleFont(((TObjString*)tok->At(0))->GetString().Atoi());
1424 ax->SetTitleSize(((TObjString*)tok->At(1))->GetString().Atof());
1425 ax->SetTitleOffset(((TObjString*)tok->At(2))->GetString().Atof());
1426 }
1427 }
1428
1429 if (!label.IsNull()) {
1430 tok = label.Tokenize(" ");
1431 if (tok->GetEntries() == 3) {
1432 ax->SetLabelFont(((TObjString*)tok->At(0))->GetString().Atoi());
1433 ax->SetLabelSize(((TObjString*)tok->At(1))->GetString().Atof());
1434 ax->SetLabelOffset(((TObjString*)tok->At(2))->GetString().Atof());
1435 }
1436 }
1437
1438 if (tok) delete tok;
1439
1440}
1441
1442
1443//###############################################################################################################"
1444//-------------------------------------------------
1445
1448
1450{
1451//-------------------------------------------------
1452
1453 TObjArray* tok = NULL;
1454 if (ob->IsA()->InheritsFrom("TAttLine")) {
1455 if (!line.IsNull()) {
1456 tok = line.Tokenize(" ");
1457 if (tok->GetEntries() == 3) {
1458 ob->SetLineColor(((TObjString*)tok->At(0))->GetString().Atoi());
1459 ob->SetLineStyle(((TObjString*)tok->At(1))->GetString().Atoi());
1460 ob->SetLineWidth(((TObjString*)tok->At(2))->GetString().Atoi());
1461 }
1462 }
1463 }
1464 if (tok) delete tok;
1465
1466}
1467
1468
1469//###############################################################################################################"
1470//-------------------------------------------------
1471
1474
1476{
1477//-------------------------------------------------
1478
1479 TObjArray* tok = NULL;
1480 if (ob->IsA()->InheritsFrom("TAttMarker")) {
1481 if (!marker.IsNull()) {
1482 tok = marker.Tokenize(" ");
1483 if (tok->GetEntries() == 3) {
1484 ob->SetMarkerColor(((TObjString*)tok->At(0))->GetString().Atoi());
1485 ob->SetMarkerStyle(((TObjString*)tok->At(1))->GetString().Atoi());
1486 ob->SetMarkerSize(((TObjString*)tok->At(2))->GetString().Atof());
1487 }
1488 }
1489 }
1490 if (tok) delete tok;
1491
1492}
1493
1494
1495//###############################################################################################################"
1496//-------------------------------------------------
1497
1500
1502{
1503//-------------------------------------------------
1504
1506 DefineMarkerStyle((TAttMarker*)ob, marker);
1507
1508}
1509
1510
1511//###############################################################################################################"
1512//-------------------------------------------------
1513
1516
1518{
1519//-------------------------------------------------
1520
1521 ob->GetXaxis()->SetTitle(xtit);
1522 ob->GetYaxis()->SetTitle(ytit);
1523
1524}
1525
1526//###############################################################################################################"
1527//-------------------------------------------------
1528
1531
1533{
1534//-------------------------------------------------
1535
1536 ob->GetXaxis()->SetTitle(xtit);
1537 ob->GetYaxis()->SetTitle(ytit);
1538
1539}
1540
1541//###############################################################################################################"
1542//-------------------------------------------------
1543
1546
1548{
1549//-------------------------------------------------
1550
1551 ob->GetXaxis()->SetTitle(xtit);
1552 ob->GetYaxis()->SetTitle(ytit);
1553
1554}
1555
1556//###############################################################################################################"
1557//-------------------------------------------------
1558
1568
1570{
1571 // Return value of abscissa X for which the interpolated value
1572 // of the histogram contents is equal to the given value, val.
1573 // We use the false position method (which should always converge...)
1574 // eps is required precision, i.e. convergence condition is that no further change
1575 // in result greater than eps is found.
1576 // nmax = maximum number of iterations
1577 // A solution is searched for X between the limits Xmin and Xmax of the X axis of the histo
1578 // unless arguments (xmin,xmax) are given to bracket the search
1579
1580 TSpline5* spline = new TSpline5(ob);
1581 Double_t Xmax;
1582 Double_t Xmin;
1583 if (xmin == xmax) {
1584 Xmax = ob->GetXaxis()->GetXmax();
1585 Xmin = ob->GetXaxis()->GetXmin();
1586 }
1587 else {
1588 Xmax = xmax;
1589 Xmin = xmin;
1590 }
1591
1592 Int_t n, side = 0;
1593 Double_t r = 0, fr, fs, s, ft, t;
1594 s = Xmin;
1595 t = Xmax;
1596 fs = spline->Eval(s) - val;
1597 ft = spline->Eval(t) - val;
1598
1599 for (n = 1; n <= nmax; n++) {
1600 r = (fs * t - ft * s) / (fs - ft);
1601 if (fabs(t - s) < eps * fabs(t + s)) break;
1602 fr = spline->Eval(r) - val;
1603
1604 if (fr * ft > 0) {
1605 t = r;
1606 ft = fr;
1607 if (side == -1) fs /= 2;
1608 side = -1;
1609 }
1610 else if (fs * fr > 0) {
1611 s = r;
1612 fs = fr;
1613 if (side == +1) ft /= 2;
1614 side = +1;
1615 }
1616 else break;
1617 }
1618 delete spline;
1619 return r;
1620}
1621
1622//###############################################################################################################"
1623//-------------------------------------------------
1624
1662
1663TF1* KVHistoManipulator::RescaleX(TH1* hist1, TH1* hist2, Int_t degree, Double_t* params,
1664 Int_t npoints, const Char_t* direction, Double_t xmin, Double_t xmax, Double_t qmin, Double_t qmax, Double_t eps)
1665{
1666 // Find the polynomial transformation of the X-axis of 1-D histogram hist1
1667 // so that the distribution ressembles that of histogram hist2.
1668 // Returns a TF1 which can be used to rescale hist1 (hist1 and hist2 are not modified).
1669 // User's responsibility to delete the TF1.
1670 // After fitting, the array params is filled with:
1671 // params[0] = a0
1672 // params[1] = a1
1673 // ...
1674 // params[degree] = an
1675 // params[degree+1] = Chisquare/NDF
1676 // Therefore params MUST BE at least of dimension (degree+2)
1677 //
1678 // npoints : by default (npoints=-1), we use npoints=degree+2 values of comparison
1679 // of the two cumulative distributions (see below).
1680 // more comparison points can be used by setting npoints>=degree+2.
1681 // direction : by default ("C") we use the cumulative histogram summed from low x to high x.
1682 // if direction="D", we use the cumulative histogram summed from high x to low x.
1683 // xmin, xmax : range of values of abscissa used to build cumulative histograms. default
1684 // (xmin=xmax=-1) include all values.
1685 // qmin, qmax : minimum & maximum values of cumulative histograms used for the
1686 // comparison (see below). by default qmin=0.05, qmax=0.95.
1687 //
1688 // METHOD
1689 // ======
1690 // 'hist1' contains the distribution P1 of variable X1
1691 // 'hist2' contains the distribution P2 of variable X2
1692 // Supposing that we can write X2=f_n(X1), with f_n(x)=a0 + a1*x + a2*x^2 + ... + an*x^n,
1693 // what are the parameters of the polynomial for which P1(f_n(X1))=P2(X2) ?
1694 //
1695 // Consider the cumulative distributions, C1(X1) and C2(X2).
1696 // For npoints=2 we compare the 2 X1 & X2 values for which C1=C2=qmin, qmax
1697 // For npoints=3 we compare the 3 X1 & X2 values for which C1=C2=qmin, qmin+(qmax-qmin)/2, qmax
1698 // For npoints=4 we compare the 4 X1 & X2 values for which C1=C2=qmin, qmin+(qmax-qmin)/3, qmin+2*(qmax-qmin)/3, qmax
1699 // etc. etc.
1700 // In each case we fit the npoints couples (X1,X2) with f_n
1701 //
1702
1703 int i;
1704 // calculate comparison points
1705 npoints = TMath::Max(npoints, degree + 2);
1706 TString func_name = Form("pol%d", degree);
1707 TF1* fonc = new TF1("f", func_name.Data());
1708 fonc->SetName(Form("RescaleX-%s", func_name.Data()));
1709 RescaleX(hist1, hist2, fonc, npoints, direction, xmin, xmax, qmin, qmax, eps);
1710 for (i = 0; i < degree + 1; i++) {
1711 params[i] = fonc->GetParameter(i);
1712 }
1713 Double_t chisquare = fonc->GetChisquare();
1714 if (fonc->GetNDF() > 0.0) chisquare /= fonc->GetNDF();
1715 params[degree + 1] = chisquare;
1716
1717 return fonc;
1718}
1719
1720//###############################################################################################################"
1721//-------------------------------------------------
1722
1788
1790 Option_t* opt, Int_t npoints, const Char_t* direction, Double_t xmin, Double_t xmax, Double_t qmin, Double_t qmax,
1791 Double_t eps)
1792{
1793 // Uses RescaleX(TH1* hist1, TH1* hist2, Int_t degree, Double_t* params, Double_t eps)
1794 // to find a polynomial transformation of 'hist1' abscissa so that its
1795 // distribution resembles that of 'hist2', then generates a new rescaled version of 'hist1'.
1796 //
1797 // degree = degree of polynomial to use
1798 // params = array of dimension [degree+2], after rescaling it will contain
1799 // the values of fitted parameters of polynomial, plus the Chi2/NDF of the fit:
1800 // params[0] = a0
1801 // params[1] = a1
1802 // ...
1803 // params[degree] = an
1804 // params[degree+1] = Chisquare/NDF = Chisquare (NDF is always equal to 1)
1805 //
1806 // npoints : by default (npoints=-1), we use npoints=degree+2 values of comparison
1807 // of the two cumulative distributions (see method RescaleX).
1808 // more comparison points can be used by setting npoints>=degree+2.
1809 //
1810 // direction : by default ("C") we use the cumulative histogram summed from low x to high x.
1811 // if direction="D", we use the cumulative histogram summed from high x to low x.
1812 // xmin, xmax : range of values of abscissa used to build cumulative histograms. default
1813 // (xmin=xmax=-1) include all values.
1814 // qmin, qmax : minimum & maximum values of cumulative histograms used for the
1815 // comparison (see method RescaleX). by default qmin=0.05, qmax=0.95.
1816 //
1817 // eps = relative precision used to find comparison points of histos (default = 1.e-07)
1818 //
1819 // OPTIONS
1820 // =======
1821 // opt = "norm" : rescaled histogram normalised to have same integral as hist2
1822 // opt = "bins" : rescaled histogram will have same number of bins & same limits as hist2
1823 // opt = "normbins" |
1824 // opt = "binsnorm" |--> rescaled histogram can be superposed and added to hist2
1825 //
1826 // EXAMPLE OF USE
1827 // =======================
1828 // In the following example, we fill two histograms with different numbers of random
1829 // values drawn from two Gaussian distributions with different centroids and widths.
1830 // We also add to each histogram a 'pedestal' peak which is unrelated to the 'physical'
1831 // distributions, and significant enough so that it does not permit a correct scaling of
1832 // the physical part of the distribution.
1833 //
1834 // We can overcome the problem of the 'pedestal' by
1835 // using the 'inverse cumulated distribution' and excluding the channels
1836 // where the pedestals are present. The correct scaling of the physical
1837 // distributions is recovered, as shown below:
1838 /*
1839 BEGIN_MACRO(source)
1840 {
1841 TH1F* h10 = new TH1F("h10","gaussian",4096,-.5,4095.5);
1842 TF1 g10("g10","gaus(0)",0,4100);
1843 g10.SetParameters(1,1095,233);
1844 h10->FillRandom("g10",56130);
1845 h10->SetBinContent(85,13425);
1846 TH1F* h20 = new TH1F("h20","gaussian",4096,-.5,4095.5);
1847 g10.SetParameters(1,1673,487);
1848 h20->FillRandom("g10",21370);
1849 h20->SetBinContent(78,17900);
1850 KVHistoManipulator HM2;
1851 HM2.SetVisDebug(); // turn on visual debugging -> create canvas showing rescaling procedure
1852 Double_t params0[10];
1853 // make rescaling using inverse cumulative distribution, limit to x=[100,4095]
1854 TH1* sc30 =HM2.MakeHistoRescaleX(h10,h20,1,params0,"binsnorm",5,"D",100,4095);
1855 return gROOT->GetListOfCanvases()->FindObject("VDCanvas");
1856 }
1857 END_MACRO
1858 */
1859
1860 TF1* scalefunc = RescaleX(hist1, hist2, degree, params, npoints, direction, xmin, xmax, qmin, qmax, eps);
1861 TString options(opt);
1862 options.ToUpper();
1863 Bool_t norm = options.Contains("NORM");
1864 Bool_t bins = options.Contains("BINS");
1865 Int_t nx = (bins ? hist2->GetNbinsX() : -1);
1866 Double_t nxmin = (bins ? hist2->GetXaxis()->GetXmin() : -1);
1867 Double_t nxmax = (bins ? hist2->GetXaxis()->GetXmax() : -1);
1868 TH1* scalehisto = ScaleHisto(hist1, scalefunc, 0, nx, -1, nxmin, nxmax, -1.0, -1.0, "width");
1869 if (norm) scalehisto->Scale(hist2->Integral("width") / scalehisto->Integral("width"));
1870 if (kVisDebug) {
1871 fVDCanvas->cd(4);
1872 scalehisto->DrawCopy()->SetLineColor(kRed);
1873 hist2->DrawCopy("same")->SetLineColor(kBlack);
1874 gPad->SetLogy(kTRUE);
1875 gPad->Modified();
1876 gPad->Update();
1877 }
1878 delete scalefunc;
1879 return scalehisto;
1880}
1881
1882//###############################################################################################################"
1883//-------------------------------------------------
1884
1909
1910void KVHistoManipulator::RescaleX(TH1* hist1, TH1* hist2, TF1* scale_func, Int_t npoints,
1911 const Char_t* direction, Double_t xmin, Double_t xmax, Double_t qmin, Double_t qmax, Double_t eps)
1912{
1913 // Find the transformation of the X-axis of 1-D histogram hist1
1914 // so that the distribution ressembles that of histogram hist2.
1915 // The user provides a function f(x) (TF1* scale_func) which is supposed to
1916 // transform the abscissa X1 of 'hist1' in such a way that P1(f(X1)) = P2(X2).
1917 // We fit 'npoints' comparison points (see below), npoints>=2.
1918 //
1919 // direction : by default ("C") we use the cumulative histogram summed from low x to high x.
1920 // if direction="D", we use the cumulative histogram summed from high x to low x.
1921 // xmin, xmax : range of values of abscissa used to build cumulative histograms. default
1922 // (xmin=xmax=-1) include all values.
1923 // qmin, qmax : minimum & maximum values of cumulative histograms used for the
1924 // comparison (see below). by default qmin=0.05, qmax=0.95.
1925 // METHOD
1926 // ======
1927 // 'hist1' contains the distribution P1 of variable X1
1928 // 'hist2' contains the distribution P2 of variable X2
1929 // Supposing that we can write X2=f(X1), we compare the abscissa of different
1930 // points of the two cumulative distributions, C1(X1) and C2(X2).
1931 // For npoints=2 we compare the 2 X1 & X2 values for which C1=C2=qmin, qmax
1932 // For npoints=3 we compare the 3 X1 & X2 values for which C1=C2=qmin, qmin+(qmax-qmin)/2, qmax
1933 // For npoints=4 we compare the 4 X1 & X2 values for which C1=C2=qmin, qmin+(qmax-qmin)/3, qmin+2*(qmax-qmin)/3, qmax
1934 // etc. etc.
1935 // In each case we fit the 'npoints' couples (X1,X2) with the TF1 pointed to by 'scale_func'
1936
1937 if (kVisDebug) {
1938 if (!fVDCanvas) fVDCanvas = new TCanvas("VDCanvas", "KVHistoManipulator::RescaleX");
1939 gStyle->SetOptStat("");
1940 fVDCanvas->Clear();
1941 fVDCanvas->Divide(2, 2);
1942 fVDCanvas->cd(1);
1943 hist1->DrawCopy()->SetLineColor(kBlue);
1944 hist2->DrawCopy("same")->SetLineColor(kBlack);
1945 gPad->SetLogy(kTRUE);
1946 gPad->Modified();
1947 gPad->Update();
1948 }
1949 int i;
1950 npoints = TMath::Max(2, npoints);
1951 Info("RescaleX", "Calculating transformation of histo %s using reference histo %s, %d points of comparison",
1952 hist1->GetName(), hist2->GetName(), npoints);
1953 TH1* cum1 = 0;
1954 TH1* cum2 = 0;
1955 if (xmin > -1 && xmax > -1) {
1956 cum1 = CumulatedHisto(hist1, xmin, xmax, direction, "max");
1957 cum2 = CumulatedHisto(hist2, xmin, xmax, direction, "max");
1958 }
1959 else {
1960 cum1 = CumulatedHisto(hist1, direction, -1, -1, "max");
1961 cum2 = CumulatedHisto(hist2, direction, -1, -1, "max");
1962 }
1963 if (kVisDebug) {
1964 fVDCanvas->cd(2);
1965 cum1->DrawCopy()->SetLineColor(kBlue);
1966 cum2->DrawCopy("same")->SetLineColor(kBlack);
1967 gPad->Modified();
1968 gPad->Update();
1969 }
1970 // calculate comparison points
1971 Double_t* quantiles = new Double_t[npoints];
1972 Double_t delta_q = (qmax - qmin) / (1.0 * (npoints - 1));
1973 for (i = 0; i < npoints; i++) quantiles[i] = qmin + i * delta_q;
1974 // get X1 and X2 values corresponding to quantiles
1975 Double_t* X1 = new Double_t[npoints];
1976 Double_t* X2 = new Double_t[npoints];
1977 for (i = 0; i < npoints; i++) {
1978 X1[i] = GetX(cum1, quantiles[i], eps);
1979 X2[i] = GetX(cum2, quantiles[i], eps);
1980 }
1981 for (i = 0; i < npoints; i++) {
1982 printf("COMPARISON: i=%d quantile=%f X1=%f X2=%f\n",
1983 i, quantiles[i], X1[i], X2[i]);
1984 }
1985 // fill TGraph with points to fit
1986 TGraph* fitgraph = new TGraph(npoints, X1, X2);
1987 TString fitoptions = "0N";
1988 if (kVisDebug) fitoptions = "";
1989 if (fitgraph->Fit(scale_func, fitoptions.Data()) != 0) {
1990 Error("RescaleX", "Fitting with function %s failed to converge",
1991 scale_func->GetName());
1992 }
1993 if (kVisDebug) {
1994 fVDCanvas->cd(3);
1995 fitgraph->SetMarkerStyle(20);
1996 gStyle->SetOptStat(1011);
1997 fitgraph->DrawClone("ap");
1998 gPad->Modified();
1999 gPad->Update();
2000 }
2001 delete cum1;
2002 delete cum2;
2003 delete fitgraph;
2004 delete [] quantiles;
2005 delete [] X1;
2006 delete [] X2;
2007}
2008
2009//###############################################################################################################"
2010//-------------------------------------------------
2011
2035
2036TH1* KVHistoManipulator::MakeHistoRescaleX(TH1* hist1, TH1* hist2, TF1* scale_func, Int_t npoints,
2037 Option_t* opt, const Char_t* direction, Double_t xmin, Double_t xmax, Double_t qmin, Double_t qmax, Double_t eps)
2038{
2039 // Uses RescaleX(TH1* hist1, TH1* hist2, TF1* scale_func, Int_t npoints, Double_t eps)
2040 // to transform 'hist1' abscissa using TF1 'scale_func' so that its
2041 // distribution resembles that of 'hist2', then generates a new rescaled version of 'hist1'.
2042 //
2043 // npoints = number of points of comparison between the two histograms.
2044 // Make sure this is sufficient for the TF1 used in the transformation.
2045 // i.e. for a polynomial of degree 1 (a+bx), 2 points are enough,
2046 // 3 will give a meaningful Chi^2 value.
2047 // direction : by default ("C") we use the cumulative histogram summed from low x to high x.
2048 // if direction="D", we use the cumulative histogram summed from high x to low x.
2049 // xmin, xmax : range of values of abscissa used to build cumulative histograms. default
2050 // (xmin=xmax=-1) include all values.
2051 // qmin, qmax : minimum & maximum values of cumulative histograms used for the
2052 // comparison (see method RescaleX). by default qmin=0.05, qmax=0.95.
2053 // eps = relative precision used to find comparison points of histos (default = 1.e-07)
2054 //
2055 // OPTIONS
2056 // =======
2057 // opt = "norm" : rescaled histogram normalised to have same integral as hist2
2058 // opt = "bins" : rescaled histogram will have same number of bins & same limits as hist2
2059 // opt = "normbins" |
2060 // opt = "binsnorm" |--> rescaled histogram can be superposed and added to hist2
2061
2062 RescaleX(hist1, hist2, scale_func, npoints, direction, xmin, xmax, qmin, qmax, eps);
2063 TString options(opt);
2064 options.ToUpper();
2065 Bool_t norm = options.Contains("NORM");
2066 Bool_t bins = options.Contains("BINS");
2067 Int_t nx = (bins ? hist2->GetNbinsX() : -1);
2068 Double_t nxmin = (bins ? hist2->GetXaxis()->GetXmin() : -1);
2069 Double_t nxmax = (bins ? hist2->GetXaxis()->GetXmax() : -1);
2070 TH1* scalehisto = ScaleHisto(hist1, scale_func, 0, nx, -1, nxmin, nxmax, -1.0, -1.0, "width");
2071 if (norm) scalehisto->Scale(hist2->Integral("width") / scalehisto->Integral("width"));
2072 if (kVisDebug) {
2073 fVDCanvas->cd(4);
2074 scalehisto->DrawCopy()->SetLineColor(kRed);
2075 hist2->DrawCopy("same")->SetLineColor(kBlack);
2076 gPad->SetLogy(kTRUE);
2077 gPad->Modified();
2078 gPad->Update();
2079 }
2080 return scalehisto;
2081}
2082
2083//-------------------------------------------------
2084
2088
2090{
2091 // Cumule le contenu de l histo hh entre xmin et xmax et retourne l histo correspondant
2092 // Voir CumulatedHisto(TH1* ,TString ,Int_t ,Int_t , Option_t*).
2093 Int_t bmin = hh->FindBin(xmin);
2094 Int_t bmax = hh->FindBin(xmax);
2095 if (bmax > hh->GetNbinsX()) bmax = hh->GetNbinsX();
2096 return CumulatedHisto(hh, direction, bmin, bmax, norm);
2097}
2098
2099
2100//-------------------------------------------------
2101
2109
2111{
2112 //Camcul du chi2 entre un histogramme et une fonction donnée
2113 //Warning : ne prend en compte que les bins avec une stat>0
2114 //de l histogramme
2115 // norm = kTRUE (default), normalise la valeur du Chi2 au nombre de bin pris en compte
2116 // err = kTRUE (default), prend en compte les erreurs du contenu des bins dans le calcul
2117 // (si celle ci est >0)
2118 Double_t chi2 = 0;
2119 Int_t nbre = 0;
2120 if (h1->InheritsFrom("TH2")) {
2121 Double_t xx[2];
2122 for (Int_t nx = 1; nx < h1->GetNbinsX(); nx += 1)
2123 for (Int_t ny = 1; ny <= h1->GetNbinsY(); ny += 1) {
2124 Double_t hval = h1->GetBinContent(nx, ny);
2125 if (hval > 0) {
2126 if (err) {
2127 Double_t herr = h1->GetBinError(nx, ny);
2128 if (herr > 0) {
2129 nbre += 1;
2130 xx[0] = h1->GetXaxis()->GetBinCenter(nx);
2131 xx[1] = h1->GetYaxis()->GetBinCenter(ny);
2132 Double_t fval = f1->EvalPar(xx, para);
2133 chi2 += TMath::Power((hval - fval) / herr, 2.);
2134 }
2135 }
2136 else {
2137 nbre += 1;
2138 xx[0] = h1->GetXaxis()->GetBinCenter(nx);
2139 xx[1] = h1->GetYaxis()->GetBinCenter(ny);
2140 Double_t fval = f1->EvalPar(xx, para);
2141 chi2 += TMath::Power((hval - fval), 2.);
2142 }
2143 }
2144 }
2145 }
2146 else {
2147 Double_t xx[1];
2148 for (Int_t nx = 1; nx < h1->GetNbinsX(); nx += 1) {
2149 Double_t hval = h1->GetBinContent(nx);
2150 if (hval > 0) {
2151 if (err) {
2152 Double_t herr = h1->GetBinError(nx);
2153 if (herr > 0) {
2154 nbre += 1;
2155 xx[0] = h1->GetXaxis()->GetBinCenter(nx);
2156 Double_t fval = f1->EvalPar(xx, para);
2157 chi2 += TMath::Power((hval - fval) / herr, 2.);
2158 }
2159 }
2160 else {
2161 nbre += 1;
2162 xx[0] = h1->GetXaxis()->GetBinCenter(nx);
2163 Double_t fval = f1->EvalPar(xx, para);
2164 chi2 += TMath::Power((hval - fval), 2.);
2165 }
2166 }
2167 }
2168 }
2169
2170 if (nbre == 0) {
2171 printf("Warning, KVHistoManipulator::GetChisquare :\n\taucune cellule price en compte dans le calcul du Chi2 ...\n");
2172 return -1;
2173 }
2174 return (norm ? chi2 / nbre : chi2);
2175
2176
2177}
2178
2179//-------------------------------------------------
2180
2188
2190{
2191 //Calcul du chi2 entre un histogramme et une fonction donnée
2192 //Warning : ne prend en compte que les bins avec une stat>0
2193 //de l histogramme
2194 // norm = kTRUE (default), normalise la valeur du Chi2 au nombre de bin pris en compte
2195 // err = kTRUE (default), prend en compte les erreurs du contenu des bins dans le calcul
2196 // (si celle ci est >0)
2197 Double_t chi2 = 0;
2198 Int_t nbre = 0;
2199 if (h1->InheritsFrom("TH2")) {
2200 Double_t xx[2];
2201 for (Int_t nx = 1; nx < h1->GetNbinsX(); nx += 1)
2202 for (Int_t ny = 1; ny <= h1->GetNbinsY(); ny += 1) {
2203 Double_t hval = h1->GetBinContent(nx, ny);
2204 if (hval > 0) {
2205 nbre += 1;
2206 xx[0] = h1->GetXaxis()->GetBinCenter(nx);
2207 xx[1] = h1->GetYaxis()->GetBinCenter(ny);
2208 Double_t fval = f1->EvalPar(xx, para);
2209 Double_t logfval = TMath::Log(fval);
2210 chi2 += fval - 1 * hval * logfval;
2211 }
2212 }
2213 }
2214 else {
2215 Double_t xx[1];
2216 for (Int_t nx = 1; nx < h1->GetNbinsX(); nx += 1) {
2217 Double_t hval = h1->GetBinContent(nx);
2218 if (hval > 0) {
2219 nbre += 1;
2220 xx[0] = h1->GetXaxis()->GetBinCenter(nx);
2221 Double_t fval = f1->EvalPar(xx, para);
2222 Double_t logfval = TMath::Log(fval);
2223 chi2 += fval - 1 * hval * logfval;
2224 }
2225 }
2226 }
2227
2228 if (nbre == 0) {
2229 printf("Warning, KVHistoManipulator::GetChisquare :\n\taucune cellule price en compte dans le calcul du Chi2 ...\n");
2230 return -1;
2231 }
2232 return (norm ? chi2 / nbre : chi2);
2233
2234}
2235
2236
2237
2244
2246{
2247 // Create and fill a TGraph containing, for each point in G1 and G2,
2248 // the value of (y1/y2).
2249 // G1 & G2 should have the same number of points, with the same x-coordinates
2250 // i.e. x1 = x2 for all points
2251 // if any y2 value = 0, we set the corresponding point's y=0
2252
2253 Int_t npoints = G1->GetN();
2254 if (G2->GetN() != npoints) {
2255 Error("DivideGraphs", "Graphs must have same number of points");
2256 return nullptr;
2257 }
2258 // make copy of G1
2259 AUTO_NEW_CTOR(TGraph, Gdiv)(*G1);
2260 Gdiv->SetName(Form("%s_divided_by_%s", G1->GetName(), G2->GetName()));
2261 Gdiv->SetTitle(Form("%s divided by %s", G1->GetTitle(), G2->GetTitle()));
2262 Double_t* X = G1->GetX();
2263 Double_t* Y1 = G1->GetY();
2264 Double_t* Y2 = G2->GetY();
2265 for (int i = 0; i < npoints; i++) {
2266 if (Y2[i] != 0) Gdiv->SetPoint(i, X[i], Y1[i] / Y2[i]);
2267 else Gdiv->SetPoint(i, X[i], 0.);
2268 }
2269 return Gdiv;
2270}
2271
2272
2273
2275
2277{
2278 Int_t npoints = g0->GetN();
2279 if (g1->GetN() != npoints) {
2280 Error("ComputeNewGraphFrom", "Graphs must have same number of points %d != %d", npoints, g1->GetN());
2281 return nullptr;
2282 }
2283
2284 TF1 f1("func_ComputeNewGraphFrom", formula, 0, 1);
2285 if (f1.IsZombie() || f1.GetNpar() != 2) {
2286 Error("ComputeNewGraphFrom", "formula %s for the operation is not valid or has not 2 parameters", formula.Data());
2287 return nullptr;
2288 }
2289
2290 auto gfinal = new TGraph;
2291 gfinal->SetName(Form("from_%s_%s", g0->GetName(), g1->GetName()));
2292 Double_t* x0 = g0->GetX();
2293 Double_t* y0 = g0->GetY();
2294
2295 Double_t* x1 = g1->GetX();
2296 Double_t* y1 = g1->GetY();
2297
2298 for (Int_t ii = 0; ii < npoints; ii++) {
2299 f1.SetParameters(y0[ii], y1[ii]);
2300 if (x1[ii] != x0[ii])
2301 Warning("ComputeNewGraphFrom", "X values are different for the same point %d : %lf %lf", ii, x0[ii], x1[ii]);
2302 Double_t result = f1.Eval(x0[ii]);
2303 gfinal->SetPoint(ii, x0[ii], result);
2304 }
2305 return gfinal;
2306}
2307
2308
2309
2321
2323{
2324
2325 //generalization of the previous method ComputeNewGraphFrom(TGraph* g0, TGraph* g1, TString formula)
2326 //
2327 //create a new TGraph from the combination of the graph in the list
2328 //new errors are not computed
2329 //x-values are unchanged
2330 //from a given list "lgr", compute expression indicated in the "formula", with the following prescription :
2331 //the "i" position of the graph in the list will correspond to the "[i]" parameters of the formula
2332 //for each points of the graphs, the formula is computed with the y-values of the graphs and the new y-value is set to the new graph
2333 //Examples
2334 //if there is 2 graphs in the list, and the formula "[0]+[1]" is set, the method will give a new TGraph with the sum of the y-values of the two TGraphs
2335
2336 Int_t ngr = lgr->GetEntries();
2337 TF1 f1("func_ComputeNewGraphFrom", formula, 0, 1);
2338 if (f1.IsZombie()) {
2339 Error("ComputeNewGraphFrom", "wrong formula %s, check the expression", formula.Data());
2340 return nullptr;
2341 }
2342 if (f1.GetNpar() != ngr) {
2343 Error("ComputeNewGraphFrom", "number of parameters (%d) of formula %s is not the same as the number of graphics (%d) in the list", f1.GetNpar(), formula.Data(), ngr);
2344 return nullptr;
2345 }
2346
2347 auto gfinal = new TGraph;
2348 gfinal->SetName("new_graph");
2349
2350 TGraph* gr;
2351 Int_t npoints = ((TGraph*)lgr->At(0))->GetN();
2352
2353 std::vector<Double_t> par(ngr);
2354 Double_t xx;
2355 for (Int_t ii = 0; ii < npoints; ii++) {
2356 for (Int_t jj = 0; jj < ngr; jj += 1) {
2357 gr = (TGraph*)lgr->At(jj);
2358 gr->GetPoint(ii, xx, par[jj]);
2359 }
2360 f1.SetParameters(par.data());
2361
2362 Double_t result = f1.Eval(xx);
2363 gfinal->SetPoint(ii, xx, result);
2364 }
2365 return gfinal;
2366}
2367
2368
2369
2375
2376std::vector<Double_t> KVHistoManipulator::GetLimits(TGraph* G1)
2377{
2378 /*
2379 xmin -> limits[0];
2380 ymin -> limits[1];
2381 xmax -> limits[2];
2382 ymax -> limits[3];
2383 */
2384 std::vector<Double_t> limits(4);
2385 Double_t xx, yy;
2386 for (Int_t ii = 0; ii < G1->GetN(); ii += 1) {
2387 G1->GetPoint(ii, xx, yy);
2388 if (ii == 0) {
2389 limits[0] = limits[2] = xx;
2390 limits[1] = limits[3] = yy;
2391 }
2392 else {
2393 if (xx < limits[0]) limits[0] = xx;
2394 if (yy < limits[1]) limits[1] = yy;
2395 if (xx > limits[2]) limits[2] = xx;
2396 if (yy > limits[3]) limits[3] = yy;
2397 }
2398 }
2399
2400 return limits;
2401}
2402
2403
2404
2410
2412{
2413 /*
2414 xmin -> limits[0];
2415 ymin -> limits[1];
2416 xmax -> limits[2];
2417 ymax -> limits[3];
2418 */
2419 std::vector<Double_t> limits(4);
2420 Double_t xx, yy;
2421 Double_t ex, ey;
2422 for (Int_t ii = 0; ii < G1->GetN(); ii += 1) {
2423 G1->GetPoint(ii, xx, yy);
2424 ex = G1->GetErrorX(ii);
2425 ey = G1->GetErrorY(ii);
2426
2427 if (ii == 0) {
2428 limits[0] = limits[2] = xx;
2429 limits[1] = limits[3] = yy;
2430 }
2431 else {
2432 if (xx - ex < limits[0]) limits[0] = xx - ex;
2433 if (yy - ey < limits[1]) limits[1] = yy - ey;
2434 if (xx + ex > limits[2]) limits[2] = xx + ex;
2435 if (yy + ey > limits[3]) limits[3] = yy + ey;
2436 }
2437 }
2438
2439 return limits;
2440
2441}
2442
2443
2444
2450
2452{
2453
2454 /*
2455 xmin -> limits[0];
2456 ymin -> limits[1];
2457 xmax -> limits[2];
2458 ymax -> limits[3];
2459 */
2460 std::vector<Double_t> limits, temp;
2461
2462 TList* lg = mgr->GetListOfGraphs();
2463 TGraph* gr = nullptr;
2464 for (Int_t ii = 0; ii < lg->GetEntries(); ii += 1) {
2465 gr = dynamic_cast<TGraph*>(lg->At(ii));
2466 if (ii == 0) {
2467 limits = GetLimits(gr);
2468 }
2469 else {
2470 temp = GetLimits(gr);
2471 if (temp[0] < limits[0]) limits[0] = temp[0];
2472 if (temp[1] < limits[1]) limits[1] = temp[1];
2473 if (temp[2] > limits[2]) limits[2] = temp[2];
2474 if (temp[3] > limits[3]) limits[3] = temp[3];
2475 }
2476 }
2477
2478 return limits;
2479
2480}
2481
2482
2483
2489
2490std::vector<Double_t> KVHistoManipulator::GetLimits(TProfile* G1)
2491{
2492 /*
2493 xmin -> limits[0];
2494 ymin -> limits[1];
2495 xmax -> limits[2];
2496 ymax -> limits[3];
2497 */
2498 std::vector<Double_t> limits(4);
2499 Double_t xx, yy;
2500 Bool_t first = kTRUE;
2501 for (Int_t ii = 1; ii <= G1->GetNbinsX(); ii += 1) {
2502 Double_t stat = G1->GetBinEntries(ii);
2503 if (stat > 0) {
2504 xx = G1->GetBinCenter(ii);
2505 yy = G1->GetBinContent(ii);
2506 if (first) {
2507 first = kFALSE;
2508 limits[0] = limits[2] = xx;
2509 limits[1] = limits[3] = yy;
2510 }
2511 else {
2512 if (xx < limits[0]) limits[0] = xx;
2513 if (yy < limits[1]) limits[1] = yy;
2514 if (xx > limits[2]) limits[2] = xx;
2515 if (yy > limits[3]) limits[3] = yy;
2516 }
2517 }
2518 }
2519
2520 return limits;
2521
2522}
2523
2524
2525
2531
2533{
2534
2535 /*
2536 xmin -> limits[0];
2537 ymin -> limits[1];
2538 xmax -> limits[2];
2539 ymax -> limits[3];
2540 */
2541 std::vector<Double_t> temp, limits;
2542
2543 TProfile* gr = nullptr;
2544 for (Int_t ii = 0; ii < mgr->GetEntries(); ii += 1) {
2545 gr = dynamic_cast<TProfile*>(mgr->At(ii));
2546 if (ii == 0) {
2547 limits = GetLimits(gr);
2548 }
2549 else {
2550 temp = GetLimits(gr);
2551 if (temp[0] < limits[0]) limits[0] = temp[0];
2552 if (temp[1] < limits[1]) limits[1] = temp[1];
2553 if (temp[2] > limits[2]) limits[2] = temp[2];
2554 if (temp[3] > limits[3]) limits[3] = temp[3];
2555 }
2556 }
2557
2558 return limits;
2559
2560}
2561
2562
2563
2567
2569{
2570
2571 //Getthe limits of the histogram in the current pad
2572 //and apply them to the others histogram drawn on the others pads
2573
2574 if (!gPad) return;
2575 TObject* obj = nullptr;
2576 TVirtualPad* tmp = gPad;
2577 TIter nextp(gPad->GetListOfPrimitives());
2578 TH1* h1 = nullptr;
2579 while ((obj = nextp())) {
2580 if (obj->InheritsFrom("TH1")) {
2581 h1 = dynamic_cast<TH1*>(obj);
2582 }
2583 }
2584 if (h1) {
2585 Int_t x1 = h1->GetXaxis()->GetFirst();
2586 Int_t x2 = h1->GetXaxis()->GetLast();
2587 Double_t y1, y2;
2588 Double_t z1, z2;
2589
2590 if (h1->GetDimension() == 1) { //TH1 ou TProfile
2591 y1 = h1->GetMinimum();
2592 y2 = h1->GetMaximum();
2593 }
2594 else {
2595 y1 = h1->GetYaxis()->GetFirst();
2596 y2 = h1->GetYaxis()->GetLast();
2597 }
2598
2599 if (h1->GetDimension() == 2) { //TH2 ou TProfile2D
2600 z1 = h1->GetMinimum();
2601 z2 = h1->GetMaximum();
2602 }
2603 else {
2604 z1 = h1->GetZaxis()->GetFirst();
2605 z2 = h1->GetZaxis()->GetLast();
2606 }
2607
2608 std::cout << x1 << " " << x2 << " - " << y1 << " " << y2 << " - " << z1 << " " << z2 << std::endl;
2609
2610 Int_t nc = 1;
2611 TCanvas* cc = gPad->GetCanvas();
2612 TVirtualPad* pad = nullptr;
2613 while ((pad = cc->GetPad(nc))) {
2614 if (tmp != pad) {
2615 pad->cd();
2616 if (AlsoLog) {
2617 gPad->SetLogx(tmp->GetLogx());
2618 gPad->SetLogy(tmp->GetLogy());
2619 gPad->SetLogz(tmp->GetLogz());
2620 }
2621 TIter nextq(gPad->GetListOfPrimitives());
2622 TH1* h1 = nullptr;
2623 while ((obj = nextq())) {
2624 if (obj->InheritsFrom("TH1")) {
2625 h1 = dynamic_cast<TH1*>(obj);
2626
2627 h1->GetXaxis()->SetRange(x1, x2);
2628 if (h1->GetDimension() == 1) {
2629 h1->SetMinimum(y1);
2630 h1->SetMaximum(y2);
2631 }
2632 else {
2633 h1->GetYaxis()->SetRange(y1, y2);
2634 }
2635 if (h1->GetDimension() == 2) {
2636 h1->SetMinimum(z1);
2637 h1->SetMaximum(z2);
2638 }
2639 else {
2640 h1->GetZaxis()->SetRange(y1, y2);
2641 }
2642 }
2643 }
2644 gPad->Update();
2645 }
2646 nc += 1;
2647 }
2648 cc->Paint();
2649 cc->Update();
2650 }
2651 gPad = tmp;
2652}
2653
2654
int Int_t
ROOT::R::TRInterface & r
TObject * clone(const char *newname) const override
bool Bool_t
char Char_t
constexpr Bool_t kFALSE
double Double_t
constexpr Bool_t kTRUE
const char Option_t
kRed
kBlack
kBlue
#define X(type, name)
#define gDirectory
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char mode
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t width
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize fs
Option_t Option_t TPoint TPoint const char y1
char name[80]
float xmin
float ymin
float xmax
float ymax
R__EXTERN TRandom * gRandom
char * Form(const char *fmt,...)
R__EXTERN TStyle * gStyle
#define gPad
Toolkit for various operations on histograms & graphs not provided by ROOT.
TF1 * RescaleX(TH1 *hist1, TH1 *hist2, Int_t degree, Double_t *params, Int_t npoints=-1, const Char_t *direction="C", Double_t xmin=-1, Double_t xmax=-1, Double_t qmin=0.05, Double_t qmax=0.95, Double_t eps=1.e-07)
KVHistoManipulator()
Default constructor.
void DefineMarkerStyle(TAttMarker *ob, TString marker)
TH1 * MakeHistoRescaleX(TH1 *hist1, TH1 *hist2, Int_t degree, Double_t *params, Option_t *opt="", Int_t npoints=-1, const Char_t *direction="C", Double_t xmin=-1, Double_t xmax=-1, Double_t qmin=0.05, Double_t qmax=0.95, Double_t eps=1.e-07)
Bool_t kVisDebug
= kTRUE for visual debugging
TH2 * CentreeReduiteX(TH2 *hh, Int_t nx=-1, Double_t xmin=-1., Double_t xmax=-1.)
std::vector< Double_t > GetLimits(TGraph *G1)
TGraph * ComputeNewGraphFrom(TGraph *g0, TGraph *g1, const TString &formula)
TH2 * CentreeReduiteY(TH2 *hh, Int_t ny=-1, Double_t ymin=-1., Double_t ymax=-1.)
void DefineStyle(TObject *ob, TString line, TString marker)
Double_t GetX(TH1 *ob, Double_t val, Double_t eps=1.e-07, Int_t nmax=50, Double_t xmin=-1.0, Double_t xmax=-1.0)
Double_t GetChisquare(TH1 *h1, TF1 *f1, Bool_t norm=kTRUE, Bool_t err=kTRUE, Double_t *para=nullptr)
TGraph * ScaleGraph(const TGraph *hh, TF1 *fx=nullptr, TF1 *fy=nullptr) const
KVNumberList * Saucisson(TH1 *hh, Int_t ntranches=10)
TGraph * LinkGraphs(TGraph *grx, TGraph *gry)
void DefinePattern(TH1 *ob, TString titleX="42 0.08 0.8", TString titleY="42 0.07 1.2", TString labelX="42 0.05 0.005", TString labelY="42 0.05 0.006")
virtual ~KVHistoManipulator(void)
TH2 * RenormaliseHisto(TH2 *hh, Int_t bmin=-1, Int_t bmax=-1, TString axis="X", Double_t valref=1)
Int_t CutStatBin(TH1 *hh, Int_t stat_min=-1, Int_t stat_max=-1)
TGraphErrors * GetMomentEvolution(TH2 *hh, TString momentx, TString momenty, TString axis="Y", Double_t stat_min=0)
void DefineTitle(TH1 *ob, TString xtit, TString ytit)
TGraph * ExtractMeanAndSigmaFromProfile(TProfile *pf, TGraph *&sigma)
void ApplyCurrentLimitsToAllCanvas(Bool_t AlsoLog=kFALSE)
TH1 * CentreeReduite(TH1 *hh, Int_t nx=-1, Int_t ny=-1, Double_t xmin=-1., Double_t xmax=-1., Double_t ymin=-1., Double_t ymax=-1.)
TCanvas * fVDCanvas
used for visual debugging
TH1 * ScaleHisto(TH1 *hh, TF1 *fx, TF1 *fy=NULL, Int_t nx=-1, Int_t ny=-1, Double_t xmin=-1., Double_t xmax=-1., Double_t ymin=-1., Double_t ymax=-1., Option_t *norm="")
TGraph * DivideGraphs(TGraph *G1, TGraph *G2)
Double_t GetLikelihood(TH1 *h1, TF1 *f1, Bool_t norm=kTRUE, Double_t *para=nullptr)
KVList * Give_ProjectionList(TH2 *hh, Double_t MinIntegral=-1, TString axis="X")
TGraphErrors * MakeGraphFrom(TProfile *pf, Bool_t Error=kTRUE)
TH1 * CumulatedHisto(TH1 *hh, TString direction="C", Int_t bmin=-1, Int_t bmax=-1, Option_t *norm="surf")
TH1 * GetDerivative(TH1 *hh, Int_t order)
void DefineLineStyle(TAttLine *ob, TString line)
Double_t GetCorrelationFactor(TH2 *hh)
Int_t Apply_TCutG(TH2 *hh, TCutG *cut, TString mode="in")
Extended TList class which owns its objects by default.
Definition KVList.h:28
Strings used to represent a set of ranges of values.
Bool_t End(void) const
void Begin(void) const
void Add(Int_t)
Add value 'n' to the list.
Int_t Next(void) const
virtual void Add(TObject *obj)
virtual void SetTitleOffset(Float_t offset=1)
virtual void SetLabelSize(Float_t size=0.04)
virtual void SetTitleFont(Style_t font=62)
virtual void SetLabelOffset(Float_t offset=0.005)
virtual void SetLabelFont(Style_t font=62)
virtual void SetTitleSize(Float_t size=0.04)
virtual void SetLineStyle(Style_t lstyle)
virtual void SetLineWidth(Width_t lwidth)
virtual void SetLineColor(Color_t lcolor)
virtual TClass * IsA() const
virtual void SetMarkerColor(Color_t mcolor=1)
virtual void SetMarkerStyle(Style_t mstyle=1)
virtual void SetMarkerSize(Size_t msize=1)
virtual TClass * IsA() const
virtual Int_t FindBin(const char *label)
virtual Double_t GetBinCenter(Int_t bin) const
Double_t GetXmax() const
virtual Double_t GetBinLowEdge(Int_t bin) const
Int_t GetLast() const
Double_t GetXmin() const
Int_t GetNbins() const
virtual void SetRange(Int_t first=0, Int_t last=0)
virtual Double_t GetBinWidth(Int_t bin) const
virtual Double_t GetBinUpEdge(Int_t bin) const
Int_t GetFirst() const
void Clear(Option_t *option="") override
TVirtualPad * cd(Int_t subpadnumber=0) override
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
Bool_t InheritsFrom(const char *cl) const override
virtual Int_t GetEntries() const
virtual Double_t Derivative(Double_t x, Double_t *params=nullptr, Double_t epsilon=0.001) const
virtual Int_t GetNDF() const
virtual Double_t GetParameter(const TString &name) const
TAxis * GetYaxis() const
Double_t GetChisquare() const
virtual void SetRange(Double_t xmin, Double_t xmax)
virtual Int_t GetNpar() const
virtual void SetNpx(Int_t npx=100)
virtual Double_t EvalPar(const Double_t *x, const Double_t *params=nullptr)
virtual void SetParameters(const Double_t *params)
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
TAxis * GetXaxis() const
Double_t GetErrorY(Int_t bin) const override
Double_t GetErrorX(Int_t bin) const override
virtual void SetPointError(Double_t ex, Double_t ey)
virtual Int_t IsInside(Double_t x, Double_t y) const
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Double_t * GetY() const
Int_t GetN() const
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Axis_t xmin=0, Axis_t xmax=0)
Double_t * GetX() const
void SetName(const char *name="") override
TAxis * GetXaxis() const
TAxis * GetYaxis() const
virtual Double_t * GetEY() const
void SetNameTitle(const char *name="", const char *title="") override
virtual Int_t GetPoint(Int_t i, Double_t &x, Double_t &y) const
virtual Double_t GetBinCenter(Int_t bin) const
TAxis * GetZaxis()
void SetTitle(const char *title) override
virtual Int_t GetNbinsY() const
virtual Double_t GetBinError(Int_t bin) const
virtual Double_t GetMean(Int_t axis=1) const
virtual Int_t GetDimension() const
void SetNameTitle(const char *name, const char *title) override
TAxis * GetXaxis()
virtual Double_t GetMaximum(Double_t maxval=FLT_MAX) const
virtual Int_t GetNbinsX() const
virtual void SetMaximum(Double_t maximum=-1111)
virtual void SetBinError(Int_t bin, Double_t error)
TAxis * GetYaxis()
virtual Int_t Fill(const char *name, Double_t w)
virtual void SetMinimum(Double_t minimum=-1111)
Double_t GetRMS(Int_t axis=1) const
virtual void SetBinContent(Int_t bin, Double_t content)
virtual Double_t GetBinLowEdge(Int_t bin) const
void SetName(const char *name) override
virtual TH1 * DrawCopy(Option_t *option="", const char *name_postfix="_copy") const
virtual Double_t Integral(Int_t binx1, Int_t binx2, Option_t *option="") const
virtual Double_t GetBinContent(Int_t bin) const
virtual void SetBins(Int_t nx, const Double_t *xBins)
virtual Double_t GetBinWidth(Int_t bin) const
virtual void Scale(Double_t c1=1, Option_t *option="")
virtual Int_t FindBin(Double_t x, Double_t y=0, Double_t z=0)
TObject * Clone(const char *newname="") const override
virtual Double_t GetMinimum(Double_t minval=-FLT_MAX) const
TH1D * ProjectionY(const char *name="_py", Int_t firstxbin=0, Int_t lastxbin=-1, Option_t *option="") const
void SetBinContent(Int_t bin, Double_t content) override
TH1D * ProjectionX(const char *name="_px", Int_t firstybin=0, Int_t lastybin=-1, Option_t *option="") const
virtual Double_t GetBinContent(Int_t bin) const
TObject * At(Int_t idx) const override
void InitWithPrototype(const char *function, const char *proto, ROOT::EFunctionMatchMode mode=ROOT::kConversionMatch)
Bool_t IsValid() const
void Execute()
TList * GetListOfGraphs() const
virtual void SetTitle(const char *title="")
const char * GetName() const override
const char * GetTitle() const override
virtual void SetName(const char *name)
Int_t GetEntries() const override
TObject * At(Int_t idx) const override
virtual const char * ClassName() const
virtual TObject * DrawClone(Option_t *option="") const
R__ALWAYS_INLINE Bool_t IsZombie() const
virtual Bool_t InheritsFrom(const char *classname) const
void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0) override
Double_t GetBinContent(Int_t bin) const override
virtual void SetBinEntries(Int_t bin, Double_t w)
virtual Double_t GetBinEntries(Int_t bin) const
Double_t GetBinError(Int_t bin) const override
virtual Double_t Uniform(Double_t x1, Double_t x2)
virtual TObject * At(Int_t idx) const=0
Double_t Eval(Double_t x) const override
const char * Data() const
void ToUpper()
TObjArray * Tokenize(const TString &delim) const
Bool_t IsNull() const
void Form(const char *fmt,...)
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
void SetOptStat(Int_t stat=1)
virtual Int_t GetLogz() const=0
virtual TVirtualPad * cd(Int_t subpadnumber=0)=0
virtual Int_t GetLogy() const=0
virtual Int_t GetLogx() const=0
TLine * line
RVec< PromoteType< T > > abs(const RVec< T > &v)
RVec< PromoteTypes< T0, T1 > > pow(const RVec< T0 > &v, const T1 &y)
const Double_t sigma
Double_t ey[n]
const Int_t n
TGraphErrors * gr
Double_t ex[n]
TH1F * h1
TF1 * f1
TH1 * h
Expr< UnaryOp< Sqrt< T >, Expr< A, T, D, D2, R >, T >, T, D, D2, R > sqrt(const Expr< A, T, D, D2, R > &rhs)
Expr< UnaryOp< Fabs< T >, Expr< A, T, D, D2, R >, T >, T, D, D2, R > fabs(const Expr< A, T, D, D2, R > &rhs)
void Error(const char *location, const char *fmt,...)
void Info(const char *location, const char *fmt,...)
void Warning(const char *location, const char *fmt,...)
Double_t Power(Double_t x, Double_t y)
Double_t Log(Double_t x)
Double_t Sqrt(Double_t x)
Double_t Abs(Double_t d)
Double_t Max(Double_t a, Double_t b)
ClassImp(TPyArg)