CbmRoot
Loading...
Searching...
No Matches
PairAnalysisHelper.cxx
Go to the documentation of this file.
1
2// //
3// helper functions wrapped in a namespace.
4//
5//
6// Authors:
7// * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
8// Julian Book <Julian.Book@cern.ch>
9//
10// //
12
13
14#include "PairAnalysisHelper.h"
15
16#include "CbmModuleList.h"
17
18#include <TDatabasePDG.h>
19#include <TError.h>
20#include <TF1.h>
21#include <TFormula.h>
22#include <TGraph.h>
23#include <TH1.h>
24#include <TH2.h>
25#include <THashList.h>
26#include <TMCProcess.h>
27#include <TMath.h>
28#include <TObjArray.h>
29#include <TObjString.h>
30#include <TPad.h>
31#include <TParticlePDG.h>
32#include <TProfile.h>
33#include <TRandom.h>
34#include <TVectorD.h>
35
36#include "PairAnalysisStyler.h"
38
39//_____________________________________________________________________________
40TVectorD* PairAnalysisHelper::MakeLogBinning(Int_t nbinsX, Double_t xmin, Double_t xmax)
41{
42 //
43 // Make logarithmic binning
44 // the user has to delete the array afterwards!!!
45 //
46
47 //check limits
48 if (xmin < 1e-20 || xmax < 1e-20) {
49 Error("PairAnalysisHelper::MakeLogBinning", "For Log binning xmin and xmax must be > 1e-20. Using linear binning "
50 "instead!");
51 return PairAnalysisHelper::MakeLinBinning(nbinsX, xmin, xmax);
52 }
53 if (xmax < xmin) {
54 Double_t tmp = xmin;
55 xmin = xmax;
56 xmax = tmp;
57 }
58 TVectorD* binLim = new TVectorD(nbinsX + 1);
59 Double_t first = xmin;
60 Double_t last = xmax;
61 Double_t expMax = TMath::Log(last / first);
62 for (Int_t i = 0; i < nbinsX + 1; ++i) {
63 (*binLim)[i] = first * TMath::Exp(expMax / nbinsX * (Double_t) i);
64 }
65 return binLim;
66}
67
68//_____________________________________________________________________________
69TVectorD* PairAnalysisHelper::MakeLinBinning(Int_t nbinsX, Double_t xmin, Double_t xmax)
70{
71 //
72 // Make linear binning
73 // the user has to delete the array afterwards!!!
74 //
75 if (xmax < xmin) {
76 Double_t tmp = xmin;
77 xmin = xmax;
78 xmax = tmp;
79 }
80 TVectorD* binLim = new TVectorD(nbinsX + 1);
81 Double_t first = xmin;
82 Double_t last = xmax;
83 Double_t binWidth = (last - first) / nbinsX;
84 for (Int_t i = 0; i < nbinsX + 1; ++i) {
85 (*binLim)[i] = first + binWidth * (Double_t) i;
86 }
87 return binLim;
88}
89
90//_____________________________________________________________________________
92{
93 //
94 // Make arbitrary binning, bins separated by a ','
95 //
96 TString limits(bins);
97 if (limits.IsNull()) {
98 Error("PairAnalysisHelper::MakeArbitraryBinning", "Bin Limit string is empty, cannot add the variable");
99 return 0x0;
100 }
101
102 TObjArray* arr = limits.Tokenize(",");
103 Int_t nLimits = arr->GetEntries();
104 if (nLimits < 2) {
105 Error("PairAnalysisHelper::MakeArbitraryBinning", "Need at leas 2 bin limits, cannot add the variable");
106 delete arr;
107 return 0x0;
108 }
109
110 TVectorD* binLimits = new TVectorD(nLimits);
111 for (Int_t iLim = 0; iLim < nLimits; ++iLim) {
112 (*binLimits)[iLim] = (static_cast<TObjString*>(arr->At(iLim)))->GetString().Atof();
113 }
114
115 delete arr;
116 return binLimits;
117}
118
119//_____________________________________________________________________________
120TVectorD* PairAnalysisHelper::MakeGausBinning(Int_t nbinsX, Double_t mean, Double_t sigma)
121{
122 //
123 // Make gaussian binning
124 // the user has to delete the array afterwards!!!
125 //
126
127 TVectorD* binLim = new TVectorD(nbinsX + 1);
128
129 TF1 g("g", "gaus", mean - 5 * sigma, mean + 5 * sigma);
130 g.SetParameters(1, mean, sigma);
131 Double_t sum = g.Integral(mean - 5 * sigma, mean + 5 * sigma);
132 //printf("full integral gaussian: %f \n",sum);
133
134 TF1 g2("g2", "gaus(0)/[3]", mean - 5 * sigma, mean + 5 * sigma);
135 g2.SetParameters(1, mean, sigma, sum);
136
137 // Double_t *params=g2.GetParameters();
138
139 Double_t epsilon = sigma / 10000; // step size
140 Double_t xt = mean - 5 * sigma; // bin limit
141 Double_t pint = 0.0; // previous integral
142
143 Int_t bin = 0; // current entry
144 Double_t lim = epsilon; // requested integral values (start,..., end values)
145
146 // calculate intergral until you found all limits
147 while ((xt += epsilon) <= mean + 5 * sigma) {
148
149 // current integral
150 // Double_t cint = g2.Integral(mean-5*sigma,xt,params,epsilon); //fast, but NOT root 6 (Integral(min,max,epsilon)) compatible
151 Double_t cint = g2.Integral(mean - 5 * sigma, xt); //,params,epsilon); //slow
152 // printf(" integral to %f: %f , search limit: %f \n",xt,cint,lim);
153
155 if (cint >= lim && pint < lim) {
156 // printf(" %d bin found for %f with requested integral %f, actual integral is %f \n",bin,xt,lim,cint);
158 (*binLim)[bin] = xt;
159
161 bin++;
162 lim = bin * (1. / nbinsX);
164 if (bin == nbinsX) lim = 1. - epsilon;
165 }
166
168 pint = cint;
169 }
170
171 // binLim->Print();
172 return binLim;
173}
174
175//_____________________________________________________________________________
176TVectorD* PairAnalysisHelper::CombineBinning(TVectorD* low, TVectorD* high)
177{
178 //
179 // Make a new combined binning of "low" and "high"
180 // the user has to delete the returned array afterwards!!!
181 //
182
183 // fill final vector
184 Int_t lastEqualsFirst = (Int_t)(TMath::Abs((*low)[low->GetNrows() - 1] - (*high)[0]) < 1.e-15);
185 TVectorD* binLim = new TVectorD(low->GetNrows() + high->GetNrows() - lastEqualsFirst);
186 for (Int_t i = 0; i < low->GetNrows(); i++)
187 (*binLim)[i] = (*low)[i];
188 for (Int_t i = 0; i < high->GetNrows() - lastEqualsFirst; i++)
189 (*binLim)[low->GetNrows() + i] = (*high)[i + lastEqualsFirst];
190
191 // clean up
192 // delete low; delete high;
193 // binLim->Print();
194 return binLim;
195}
196
197//_____________________________________________________________________________
198TArrayD* PairAnalysisHelper::MakeStatBinLimits(TH1* h, Double_t stat)
199{
200 //
201 // get bin limits for stat. error less than 'stat'
202 //
203 if (!h || stat > 1.) return 0x0;
204
205 Double_t cont = 0.0;
206 Double_t err = 0.0;
207 Double_t from = h->GetBinLowEdge(1);
208 Double_t to = 0.0;
209
210 TArrayD* vBins = new TArrayD(1 + 1);
211 vBins->AddAt(from, 0);
212 Int_t vEle = 1;
213
214 for (Int_t i = 1; i <= h->GetNbinsX(); i++) {
215
216 if (h->GetBinContent(i) == 0.0 && h->GetBinError(i) == 0.) continue;
217
218 to = h->GetBinLowEdge(i + 1);
219 cont += h->GetBinContent(i);
220 err += (h->GetBinError(i) * h->GetBinError(i));
221 vBins->AddAt(to, vEle);
222
223 // Printf("cont %f/%f=%f err %f(%f) sum of -> rel err %f%% (current: %f%%)",
224 // h->GetBinContent(i),h->Integral(),h->GetBinContent(i)/h->Integral(), h->GetBinError(i), TMath::Sqrt(h->GetBinContent(i)),h->GetBinError(i)/h->GetBinContent(i)*100, TMath::Sqrt(err)/cont*100);
225
226 // check for new bin
227 if (TMath::Sqrt(err) / cont <= stat
228 // || (h->GetBinContent(i)/h->Integral()) < 0.01
229 // || (h->GetBinContent(i)==0.0 && h->GetBinError(i)==0. && vEle==1)
230 ) {
231 // Printf("bin from %f to %f with err %f%%",from,to,TMath::Sqrt(err)/cont*100);
232 err = 0.0;
233 cont = 0.0;
234 vEle++;
235 from = to;
236 vBins->Set(vEle + 1);
237 }
238 }
239
240 vBins->AddAt(h->GetXaxis()->GetXmax(), vBins->GetSize() - 1);
241
242 // for(Int_t i=0;i<vBins->GetSize();i++) Printf("%d %f",i,vBins->At(i));
243 return vBins;
244}
245
246//_____________________________________________________________________________
248{
249 //
250 // Make arbitrary binning using defined PDG codes
251 //
252
253 // array of pdgcodes stored in TDatabasePDG
254 TDatabasePDG* pdg = new TDatabasePDG(); //::Instance();
255 pdg->ReadPDGTable();
256 TIter next(pdg->ParticleList());
257 TGraph gr;
258 TParticlePDG* p;
259 Int_t i = 0;
260 while ((p = (TParticlePDG*) next())) {
261 if (TMath::Abs(p->PdgCode()) < 1e+6) {
262 // printf("%s -> %d \n",p->GetName(),p->PdgCode());
263 gr.SetPoint(i++, p->PdgCode(), 1.);
264 }
265 }
266 gr.Sort();
267 TVectorD* vec = new TVectorD(gr.GetN(), gr.GetX());
268 // vec->Print();
269 delete pdg;
270 return vec;
271}
272
273//_____________________________________________________________________________
274Double_t PairAnalysisHelper::EvalFormula(TFormula* form, const Double_t* values)
275{
276 //
277 // evaluate return value for formula with current parameter values
278 //
279 Double_t params[10] = {0.};
280 //put parameter values into array using variables stored as the parameters
281 for (Int_t ip = 0; ip < form->GetNpar(); ip++)
282 params[ip] = values[(UInt_t) form->GetParameter(ip)];
283 return (form->EvalPar(0x0, params));
284}
285
286//_____________________________________________________________________________
288{
289 //
290 // evaluate the formula in a readable way (for axis etc)
291 //
292 // TODO: add units, switch to TMathText, get ride of obsolete parentheses
293 TString tform(form->GetExpFormula());
294 // TMathText
295 // tform.ReplaceAll("*","\\cdot"); // multiplication sign
296 // tform.ReplaceAll("TMath::","\\"); // get ride of TMath::
297 // TLatex
298 tform.ReplaceAll("*", "#upoint"); // multiplication sign
299 // tform.ReplaceAll("Sqrt","sqrt"); // sqrt sign
300 tform.ReplaceAll("TMath::", ""); // get ride of TMath::
301 tform.ReplaceAll("()", ""); // remove function parenthesis
302 tform.ToLower(); // lower cases (e.g. Cos -> cos)
303 // replace parameter variables with proper labels
304 for (Int_t j = 0; j < form->GetNpar(); j++)
305 tform.ReplaceAll(Form("[%d]", j), PairAnalysisVarManager::GetValueLabel((UInt_t) form->GetParameter(j)));
306 return (tform);
307}
308
309//_____________________________________________________________________________
311{
312 //
313 // build formula key with parameter names
314 //
315 TString tform("f(");
316 for (Int_t j = 0; j < form->GetNpar(); j++) {
317 tform += PairAnalysisVarManager::GetValueName((UInt_t) form->GetParameter(j));
318 if (j != form->GetNpar() - 1) tform += ",";
319 }
320 tform += ")";
321 return (tform);
322}
323
324//_____________________________________________________________________________
325TFormula* PairAnalysisHelper::GetFormula(const char* name, const char* formula)
326{
327 //
328 // build a TFormula object
329 //
330 TString check(formula);
331 if (check.IsNull()) return 0x0;
332 TFormula* form = new TFormula(name, formula);
333 // compile function
334 if (form->Compile()) return 0x0;
335 //set parameter/variable identifier
336 for (Int_t i = 0; i < form->GetNpar(); i++) {
337 form->SetParName(i, PairAnalysisVarManager::GetValueName(form->GetParameter(i)));
338 // fUsedVars->SetBitNumber((Int_t)form->GetParameter(i),kTRUE);
339 }
340 return form;
341}
342
343//_____________________________________________________________________________
344void PairAnalysisHelper::SetPDGBinLabels(TH1* hist, Bool_t clean)
345{
346 //
347 // build formula key with parameter names
348 //
349
351 if (!hLbl) hLbl = hist;
352
353 TDatabasePDG* pdg = TDatabasePDG::Instance();
354 TAxis* xaxis = hLbl->GetXaxis();
355 for (Int_t i = 1; i < hist->GetNbinsX() + 1; i++) {
356 // printf("bin %d: low edge: %.0f --> %s \n",i,xaxis->GetBinLowEdge(i), pdg->GetParticle(xaxis->GetBinLowEdge(i))->GetName());
357 if (clean && !hist->GetBinContent(i)) continue;
358 Int_t pdgCode = (Int_t) xaxis->GetBinLowEdge(i);
359 TParticlePDG* p = pdg->GetParticle(pdgCode);
360 // xaxis->SetBinLabel(i,(p?p->GetName():""));
361 xaxis->SetBinLabel(i, (p ? GetPDGlabel(pdgCode) : ""));
362 }
363
364 if (hLbl->GetDimension() == 1) return;
365
366 TAxis* yaxis = hLbl->GetYaxis();
367 TString keyY = yaxis->GetName();
368 if (keyY.Contains("pdg", TString::kIgnoreCase)) {
369 for (Int_t i = 1; i < hist->GetNbinsY() + 1; i++) {
370 // printf("bin %d: low edge: %.0f --> %s \n",i,yaxis->GetBinLowEdge(i), pdg->GetParticle(yaxis->GetBinLowEdge(i))->GetName());
371 if (clean && !hist->GetBinContent(i)) continue;
372 Int_t pdgCode = (Int_t) yaxis->GetBinLowEdge(i);
373 TParticlePDG* p = pdg->GetParticle(pdgCode);
374 // yaxis->SetBinLabel(i,(p?p->GetName():""));
375 yaxis->SetBinLabel(i, (p ? GetPDGlabel(pdgCode) : ""));
376 }
377 }
378}
379
380//_____________________________________________________________________________
382{
383 //
384 // return the label in latex format corresponding to pdg code
385 //
386
387 TString name = TDatabasePDG::Instance()->GetParticle(pdg)->GetName();
388 name.ReplaceAll("dd_1_bar", "primary");
389 name.ReplaceAll("proton", "p");
390 // correct greek letters
391 /*
392 if(name.Contains("delta",TString::kIgnoreCase) ||
393 name.Contains("gamma",TString::kIgnoreCase) ||
394 name.Contains("sigma",TString::kIgnoreCase) ||
395 name.Contains("xi",TString::kIgnoreCase) ||
396 name.Contains("lambda",TString::kIgnoreCase) ||
397 name.Contains("omega",TString::kIgnoreCase) ||
398 name.Contains("eta",TString::kIgnoreCase) ||
399 name.Contains("tau",TString::kIgnoreCase) ||
400 name.Contains("phi",TString::kIgnoreCase) ||
401 name.Contains("eta",TString::kIgnoreCase) ||
402 name.Contains("upsilon",TString::kIgnoreCase) ||
403 name.Contains("nu",TString::kIgnoreCase) ||
404 name.Contains("pi",TString::kIgnoreCase) ||
405 name.Contains("rho",TString::kIgnoreCase) ) name.Prepend("#");
406 */
407 name.ReplaceAll("delta", "#delta");
408 name.ReplaceAll("gamma", "#gamma");
409 name.ReplaceAll("psi", "#psi");
410 name.ReplaceAll("sigma", "#sigma");
411 name.ReplaceAll("xi", "#xi");
412 name.ReplaceAll("lambda", "#lambda");
413 name.ReplaceAll("omega", "#omega");
414 name.ReplaceAll("eta", "#eta");
415 name.ReplaceAll("tau", "#tau");
416 name.ReplaceAll("phi", "#phi");
417 name.ReplaceAll("upsilon", "#upsilon");
418 name.ReplaceAll("nu", "#nu");
419 name.ReplaceAll("mu", "#mu");
420 name.ReplaceAll("pi", "#pi");
421 name.ReplaceAll("rho", "#rho");
422
423 // correct anti particles
424 if (name.Contains("anti")) {
425 name.ReplaceAll("anti", "#bar{");
426 name.Append("}");
427 }
428 else if (name.Contains("_bar")) {
429 name.ReplaceAll("_bar", "}");
430 name.Prepend("#bar{");
431 }
432 // correct indices
433 name.ReplaceAll("+", "^{+}");
434 name.ReplaceAll("-", "^{-}");
435 name.ReplaceAll("0", "^{0}");
436 name.ReplaceAll("_s", "_{s}");
437 name.ReplaceAll("_c", "_{c}");
438 name.ReplaceAll("_b", "_{b}");
439 name.ReplaceAll("_1", "_{1}");
440 // specials
441 name.ReplaceAll("/psi", "/#psi");
442 // Printf(" %d = %s",pdg,name.Data());
443 return name;
444}
445
446//_____________________________________________________________________________
448{
449 //
450 // build formula key with parameter names
451 //
452 TAxis* xaxis = hist->GetXaxis();
453 for (Int_t i = 1; i < hist->GetNbinsX() + 1; i++) {
454 xaxis->SetBinLabel(i, TMCProcessName[i - 1]);
455 }
456 xaxis->LabelsOption("v"); // vertical labels
457 if (gPad) {
458 xaxis->SetTitleOffset(3.6);
459 gPad->SetTopMargin(0.01);
460 gPad->SetBottomMargin(0.4);
461 }
462}
463
464
465//_____________________________________________________________________________
467{
468 //
469 // get detector name
470 //
471 if (det == static_cast<ECbmModuleId>(9)) return ("");
472 TString name = CbmModuleList::GetModuleNameCaps(det);
473
474 // TString name = CbmModuleList::GetModuleName(det);
475 return (name);
476}
477
478//_____________________________________________________________________________
479Double_t PairAnalysisHelper::GetContentMinimum(TH1* h, Bool_t inclErr)
480{
481 //
482 // get minimum bin content of histogram (having entries)
483 //
484 Int_t bin, binx, biny, binz;
485 Int_t xfirst = h->GetXaxis()->GetFirst();
486 Int_t xlast = h->GetXaxis()->GetLast();
487 Int_t yfirst = h->GetYaxis()->GetFirst();
488 Int_t ylast = h->GetYaxis()->GetLast();
489 Int_t zfirst = h->GetZaxis()->GetFirst();
490 Int_t zlast = h->GetZaxis()->GetLast();
491 Double_t minimum = FLT_MAX, value = 0., error = 0.;
492 for (binz = zfirst; binz <= zlast; binz++) {
493 for (biny = yfirst; biny <= ylast; biny++) {
494 for (binx = xfirst; binx <= xlast; binx++) {
495 bin = h->GetBin(binx, biny, binz);
496 value = h->GetBinContent(bin);
497 error = h->GetBinError(bin);
498 // Printf(" \t hist%s bin%d value%f error%f \n",h->GetTitle(),bin,value,error);
499 if (gPad->GetLogy() && (value - error) <= 0.) continue;
500 if (error > value * 0.9) continue;
501 if (inclErr) value -= h->GetBinError(bin);
502 if (value < minimum && TMath::Abs(h->GetBinError(bin) - 1.e-15) > 1.e-15) { minimum = value; }
503 }
504 }
505 }
506 // Printf(" RETURN VALUE: hist%s %f \n",h->GetTitle(),minimum);
507 return minimum;
508}
509
510Double_t PairAnalysisHelper::GetContentMaximum(TH1* h, Bool_t inclErr)
511{
512 //
513 // get maximum bin content+error of histogram (having entries)
514 //
515 Int_t bin, binx, biny, binz;
516 Int_t xfirst = h->GetXaxis()->GetFirst();
517 Int_t xlast = h->GetXaxis()->GetLast();
518 Int_t yfirst = h->GetYaxis()->GetFirst();
519 Int_t ylast = h->GetYaxis()->GetLast();
520 Int_t zfirst = h->GetZaxis()->GetFirst();
521 Int_t zlast = h->GetZaxis()->GetLast();
522 Double_t maximum = -1. * FLT_MAX, value = 0., error = 0.;
523 for (binz = zfirst; binz <= zlast; binz++) {
524 for (biny = yfirst; biny <= ylast; biny++) {
525 for (binx = xfirst; binx <= xlast; binx++) {
526 bin = h->GetBin(binx, biny, binz);
527 value = h->GetBinContent(bin);
528 error = h->GetBinError(bin);
529 if (inclErr) value += h->GetBinError(bin);
530 if (value > maximum && TMath::Abs(error - 1.e-15) > 1.e-15) { maximum = value; }
531 }
532 }
533 }
534 return maximum;
535}
536
537Double_t PairAnalysisHelper::GetQuantile(TH1* h1, Double_t p /*=0.5*/)
538{
539 //
540 // calculates the quantile of the bin contents, p=0.5 -> Median
541 // useful functionallity for plotting 2D distibutions with some extreme outliers
542 //
543 if (p < 0.0 || p > 1.) return -1.;
544 Int_t nbinsX = h1->GetNbinsX();
545 Int_t nbinsY = h1->GetNbinsY();
546 Int_t nbinsZ = h1->GetNbinsZ();
547 Int_t nbins = (nbinsX * (nbinsY ? nbinsY : 1) * (nbinsZ ? nbinsZ : 1));
548 Int_t xbin = -1;
549 Int_t ybin = -1;
550 Int_t zbin = -1;
551 Double_t val[nbins];
552 Int_t idx[nbins];
553 Int_t nfilled = 0;
554 for (Int_t i = 1; i <= nbins; i++) {
555 h1->GetBinXYZ(i, xbin, ybin, zbin);
556 if (xbin < h1->GetXaxis()->GetFirst() || xbin > h1->GetXaxis()->GetLast()) continue;
557 if (ybin < h1->GetYaxis()->GetFirst() || ybin > h1->GetYaxis()->GetLast()) continue;
558 Double_t con = h1->GetBinContent(i);
559 Double_t err = h1->GetBinError(i);
560 if (err != 0.0) {
561 // printf("bin%d %.f+-%.f \n",i,con,err);
562 val[nfilled] = con + (h1->GetDimension() < 2 ? err : 0.0); // w or w/o err?
563 nfilled++;
564 }
565 }
566 if (nfilled == 0) return -1.;
567 TMath::Sort(nfilled, val, idx, kFALSE); // kFALSE=increasing numbers
568 Int_t pos = (Int_t)((Double_t) nfilled * p);
569 //for(int i=0; i<nfilled; i++) cout << i << " " << idx[i] << " " << val[idx[i]] << endl;
570 //printf("nfilled %d quantile %f pos %d: %f \n",nfilled,p,pos,val[idx[pos]]);
571 return val[idx[pos]];
572}
573
575{
576 //
577 // normalize slices along Y in case of a 2D histogram
578 //
579 TH2* hsum = (TH2*) h->Clone("SliceInts");
580 hsum->Reset("CE");
581 for (Int_t ix = 0; ix <= hsum->GetNbinsX() + 1; ix++) {
582 Double_t ysum = h->Integral(ix, ix, 0, hsum->GetNbinsY() + 1);
583 for (Int_t iy = 0; iy <= hsum->GetNbinsY() + 1; iy++) {
584 hsum->SetBinContent(ix, iy, ysum);
585 }
586 }
588 h->Divide(hsum);
589 delete hsum;
590}
591
592void PairAnalysisHelper::CumulateSlicesX(TH2* h, Bool_t reverse, Bool_t norm)
593{
594 //
595 // caluclate cumulative sum of bins (normalized to one)
596 // for slices along X in case of a 2D histogram
597 //
599 Int_t xstart = (reverse ? h->GetNbinsX() + 1 : 0 - 1);
600 Int_t xend = (reverse ? 0 : h->GetNbinsX() + 1 + 1);
601 Int_t xincr = (reverse ? -1 : +1);
602
603 Double_t integral = 1.;
604 for (Int_t iy = 0; iy <= h->GetNbinsY() + 1; iy++) {
605 if (norm) integral = h->Integral(0, h->GetNbinsX() + 1, iy, iy);
606
607 Double_t cumInt = 0;
608 for (Int_t ix = xstart; ix != xend; ix += xincr) {
609 cumInt += h->GetBinContent(ix, iy);
610 h->SetBinContent(ix, iy, cumInt / integral);
611 }
612 }
613}
614
615void PairAnalysisHelper::Cumulate(TH1* h, Bool_t reverse, Bool_t norm)
616{
617 //
618 // caluclate cumulative sum of bins (normalized to one)
619 //
621 Int_t xstart = (reverse ? h->GetNbinsX() + 1 : 0 - 1);
622 Int_t xend = (reverse ? 0 : h->GetNbinsX() + 1 + 1);
623 Int_t xincr = (reverse ? -1 : +1);
624
625 Double_t integral = 1.;
626 if (norm) integral = h->Integral(0, h->GetNbinsX() + 1);
627
628 Double_t cumInt = 0;
629 for (Int_t ix = xstart; ix != xend; ix += xincr) {
630 cumInt += h->GetBinContent(ix);
631 h->SetBinContent(ix, cumInt / integral);
632 }
633}
634
635
636TObject* PairAnalysisHelper::FindObjectByTitle(TObjArray* arrhist, TString ref)
637{
638 //
639 // shortcut to find a certain pair type object in array
640 //
641 for (Int_t i = 0; i < arrhist->GetEntriesFast(); i++) {
642 if (!ref.CompareTo(arrhist->UncheckedAt(i)->GetTitle())) { return arrhist->UncheckedAt(i); }
643 }
644 return 0x0;
645
646 // return ( arrhist->FindObject(Form("Pair.%s",PairAnalysis::PairClassName(type))) );
647 // TString ref=Form("Pair.%s",PairAnalysis::PairClassName(type));
648}
649
650
651//_____________________________________________________________________________
653{
654 //
655 // computes the precision of a double
656 // usefull for axis ranges etc
657
658 Bool_t bfnd = kFALSE;
659 Int_t precision = 0;
660 value *= 1e6;
661 while (!bfnd) {
662 // Printf(" value %f precision %d bfnd %d",TMath::Abs(value*TMath::Power(10,precision)), precision, bfnd);
663 bfnd = ((TMath::Abs(value * TMath::Power(10, precision)) / 1e6
664 - TMath::Floor(TMath::Abs(value * TMath::Power(10, precision)) / 1e6))
665 != 0.0
666 ? kFALSE
667 : kTRUE);
668 if (!bfnd) precision++;
669 }
670
671 // Printf("precision for %f found to be %d", value, precision);
672 return precision;
673}
ECbmModuleId
Definition CbmDefs.h:39
Char_t * gr
bool first
static TString GetModuleNameCaps(ECbmModuleId moduleId)
Data class with information on a STS local track.
static const char * GetValueName(Int_t i)
static const char * GetValueLabel(Int_t i)
Int_t GetPrecision(Double_t value)
TString GetFormulaName(TFormula *form)
TObject * FindObjectByTitle(TObjArray *arrhist, TString ref)
void SetGEANTBinLabels(TH1 *hist)
TArrayD * MakeStatBinLimits(TH1 *h, Double_t stat)
TString GetFormulaTitle(TFormula *form)
Double_t GetQuantile(TH1 *h1, Double_t p=0.5)
void CumulateSlicesX(TH2 *h, Bool_t reverse=kFALSE, Bool_t norm=kFALSE)
Double_t GetContentMinimum(TH1 *h, Bool_t inclErr=kTRUE)
Double_t GetContentMaximum(TH1 *h, Bool_t inclErr=kTRUE)
TFormula * GetFormula(const char *name, const char *formula)
TVectorD * CombineBinning(TVectorD *low, TVectorD *high)
Double_t EvalFormula(TFormula *form, const Double_t *values)
TVectorD * MakeLinBinning(Int_t nbinsX, Double_t xmin, Double_t xmax)
TVectorD * MakeGausBinning(Int_t nbinsX, Double_t mean, Double_t sigma)
void SetPDGBinLabels(TH1 *hist, Bool_t clean)
TString GetDetName(ECbmModuleId det)
void Cumulate(TH1 *h, Bool_t reverse=kFALSE, Bool_t norm=kFALSE)
TString GetPDGlabel(Int_t pdg)
TVectorD * MakeArbitraryBinning(const char *bins)
TVectorD * MakeLogBinning(Int_t nbinsX, Double_t xmin, Double_t xmax)