1 ///////////////////////////////////////////////////////////////////////////////
3 // Base class for the AliTPCCalibViewer and AliTRDCalibViewer //
4 // used for the calibration monitor //
6 // Authors: Marian Ivanov (Marian.Ivanov@cern.ch) //
7 // Jens Wiechula (Jens.Wiechula@cern.ch) //
8 // Ionut Arsene (iarsene@cern.ch) //
10 ///////////////////////////////////////////////////////////////////////////////
19 //#include <TCanvas.h>
25 #include <THashTable.h>
26 #include <TObjString.h>
27 #include <TLinearFitter.h>
31 #include <TDirectory.h>
32 #include <TFriendElement.h>
34 #include "TTreeStream.h"
35 #include "AliBaseCalibViewer.h"
37 ClassImp(AliBaseCalibViewer)
39 AliBaseCalibViewer::AliBaseCalibViewer()
43 fListOfObjectsToBeDeleted(0),
44 fTreeMustBeDeleted(0),
49 // Default constructor
53 //_____________________________________________________________________________
54 AliBaseCalibViewer::AliBaseCalibViewer(const AliBaseCalibViewer &c)
58 fListOfObjectsToBeDeleted(0),
59 fTreeMustBeDeleted(0),
64 // dummy AliBaseCalibViewer copy constructor
68 fTreeMustBeDeleted = c.fTreeMustBeDeleted;
69 fListOfObjectsToBeDeleted = c.fListOfObjectsToBeDeleted;
70 fAbbreviation = c.fAbbreviation;
71 fAppendString = c.fAppendString;
74 //_____________________________________________________________________________
75 AliBaseCalibViewer::AliBaseCalibViewer(TTree* tree)
79 fListOfObjectsToBeDeleted(0),
80 fTreeMustBeDeleted(0),
85 // Constructor that initializes the calibration viewer
88 fTreeMustBeDeleted = kFALSE;
89 fListOfObjectsToBeDeleted = new TObjArray();
91 fAppendString = ".fElements";
94 //_____________________________________________________________________________
95 AliBaseCalibViewer::AliBaseCalibViewer(const Char_t* fileName, const Char_t* treeName)
99 fListOfObjectsToBeDeleted(0),
100 fTreeMustBeDeleted(0),
106 // Constructor to initialize the calibration viewer
107 // the file 'fileName' contains the tree 'treeName'
109 fFile = new TFile(fileName, "read");
110 fTree = (TTree*) fFile->Get(treeName);
111 fTreeMustBeDeleted = kTRUE;
112 fListOfObjectsToBeDeleted = new TObjArray();
114 fAppendString = ".fElements";
117 //____________________________________________________________________________
118 AliBaseCalibViewer & AliBaseCalibViewer::operator =(const AliBaseCalibViewer & param)
121 // assignment operator - dummy
122 // not yet working!!!
125 fTreeMustBeDeleted = param.fTreeMustBeDeleted;
126 fListOfObjectsToBeDeleted = param.fListOfObjectsToBeDeleted;
127 fAbbreviation = param.fAbbreviation;
128 fAppendString = param.fAppendString;
132 //_____________________________________________________________________________
133 AliBaseCalibViewer::~AliBaseCalibViewer()
136 // AliBaseCalibViewer destructor
137 // all objects will be deleted, the file will be closed, the pictures will disappear
139 if (fTree && fTreeMustBeDeleted) {
140 fTree->SetCacheSize(0);
148 for (Int_t i = fListOfObjectsToBeDeleted->GetEntriesFast()-1; i >= 0; i--) {
149 delete fListOfObjectsToBeDeleted->At(i);
151 delete fListOfObjectsToBeDeleted;
154 //_____________________________________________________________________________
155 void AliBaseCalibViewer::Delete(Option_t* option) {
157 // Should be called from AliBaseCalibViewerGUI class only.
158 // If you use Delete() do not call the destructor.
159 // All objects (except those contained in fListOfObjectsToBeDeleted) will be deleted, the file will be closed.
162 option = option; // to avoid warnings on compiling
163 if (fTree && fTreeMustBeDeleted) {
164 fTree->SetCacheSize(0);
169 delete fListOfObjectsToBeDeleted;
172 //_____________________________________________________________________________
173 void AliBaseCalibViewer::FormatHistoLabels(TH1 *histo) const {
175 // formats title and axis labels of histo
176 // removes '.fElements'
179 TString replaceString(fAppendString.Data());
180 TString *str = new TString(histo->GetTitle());
181 str->ReplaceAll(replaceString, "");
182 histo->SetTitle(str->Data());
184 if (histo->GetXaxis()) {
185 str = new TString(histo->GetXaxis()->GetTitle());
186 str->ReplaceAll(replaceString, "");
187 histo->GetXaxis()->SetTitle(str->Data());
190 if (histo->GetYaxis()) {
191 str = new TString(histo->GetYaxis()->GetTitle());
192 str->ReplaceAll(replaceString, "");
193 histo->GetYaxis()->SetTitle(str->Data());
196 if (histo->GetZaxis()) {
197 str = new TString(histo->GetZaxis()->GetTitle());
198 str->ReplaceAll(replaceString, "");
199 histo->GetZaxis()->SetTitle(str->Data());
204 //_____________________________________________________________________________
205 TFriendElement* AliBaseCalibViewer::AddReferenceTree(const Char_t* filename, const Char_t* treename, const Char_t* refname){
207 // add a reference tree to the current tree
208 // by default the treename is 'tree' and the reference treename is 'R'
210 TFile *file = new TFile(filename);
211 fListOfObjectsToBeDeleted->Add(file);
212 TTree * tree = (TTree*)file->Get(treename);
213 return AddFriend(tree, refname);
216 //_____________________________________________________________________________
217 TString* AliBaseCalibViewer::Fit(const Char_t* drawCommand, const Char_t* formula, const Char_t* cuts,
218 Double_t & chi2, TVectorD &fitParam, TMatrixD &covMatrix){
220 // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
221 // returns chi2, fitParam and covMatrix
222 // returns TString with fitted formula
225 TString formulaStr(formula);
226 TString drawStr(drawCommand);
227 TString cutStr(cuts);
230 drawStr.ReplaceAll(fAbbreviation, fAppendString);
231 cutStr.ReplaceAll(fAbbreviation, fAppendString);
232 formulaStr.ReplaceAll(fAbbreviation, fAppendString);
234 formulaStr.ReplaceAll("++", fAbbreviation);
235 TObjArray* formulaTokens = formulaStr.Tokenize(fAbbreviation.Data());
236 Int_t dim = formulaTokens->GetEntriesFast();
238 fitParam.ResizeTo(dim);
239 covMatrix.ResizeTo(dim,dim);
241 TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
242 fitter->StoreData(kTRUE);
243 fitter->ClearPoints();
245 Int_t entries = Draw(drawStr.Data(), cutStr.Data(), "goff");
248 return new TString("An ERROR has occured during fitting!");
250 Double_t **values = new Double_t*[dim+1] ;
252 for (Int_t i = 0; i < dim + 1; i++){
254 if (i < dim) centries = fTree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff");
255 else centries = fTree->Draw(drawStr.Data(), cutStr.Data(), "goff");
257 if (entries != centries) {
260 return new TString("An ERROR has occured during fitting!");
262 values[i] = new Double_t[entries];
263 memcpy(values[i], fTree->GetV1(), entries*sizeof(Double_t));
266 // add points to the fitter
267 for (Int_t i = 0; i < entries; i++){
269 for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
270 fitter->AddPoint(x, values[dim][i], 1);
274 fitter->GetParameters(fitParam);
275 fitter->GetCovarianceMatrix(covMatrix);
276 chi2 = fitter->GetChisquare();
278 TString *preturnFormula = new TString(Form("( %e+",fitParam[0])), &returnFormula = *preturnFormula;
280 for (Int_t iparam = 0; iparam < dim; iparam++) {
281 returnFormula.Append(Form("%s*(%e)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
282 if (iparam < dim-1) returnFormula.Append("+");
284 returnFormula.Append(" )");
285 delete formulaTokens;
287 for (Int_t i = 0; i < dim + 1; i++) delete [] values[i];
289 return preturnFormula;
292 //_____________________________________________________________________________
293 Double_t AliBaseCalibViewer::GetLTM(Int_t n, Double_t *array, Double_t *sigma, Double_t fraction){
295 // returns the LTM and sigma
297 Double_t *ddata = new Double_t[n];
298 Double_t mean = 0, lsigma = 0;
300 for (UInt_t i = 0; i < (UInt_t)n; i++) {
301 ddata[nPoints]= array[nPoints];
304 Int_t hh = TMath::Min(TMath::Nint(fraction * nPoints), Int_t(n));
305 AliMathBase::EvaluateUni(nPoints, ddata, mean, lsigma, hh);
306 if (sigma) *sigma = lsigma;
311 //_____________________________________________________________________________
312 Int_t AliBaseCalibViewer::GetBin(Float_t value, Int_t nbins, Double_t binLow, Double_t binUp){
313 // Returns the 'bin' for 'value'
314 // The interval between 'binLow' and 'binUp' is divided into 'nbins' equidistant bins
315 // avoid index out of bounds error: 'if (bin < binLow) bin = binLow' and vice versa
317 GetBin(value) = #frac{nbins - 1}{binUp - binLow} #upoint (value - binLow) +1
321 Int_t bin = TMath::Nint( (Float_t)(value - binLow) / (Float_t)(binUp - binLow) * (nbins-1) ) + 1;
322 // avoid index out of bounds:
323 if (value < binLow) bin = 0;
324 if (value > binUp) bin = nbins + 1;
329 //_____________________________________________________________________________
330 TH1F* AliBaseCalibViewer::SigmaCut(TH1F *histogram, Float_t mean, Float_t sigma, Float_t sigmaMax,
331 Float_t sigmaStep, Bool_t pm) {
333 // Creates a cumulative histogram Begin_Latex S(t, #mu, #sigma) End_Latex, where you can see, how much of the data are inside sigma-intervals around the mean value
334 // The data of the distribution Begin_Latex f(x, #mu, #sigma) End_Latex are given in 'histogram'
335 // 'mean' and 'sigma' are Begin_Latex #mu End_Latex and Begin_Latex #sigma End_Latex of the distribution in 'histogram', to be specified by the user
336 // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, Begin_Latex t #sigma End_Latex)
337 // sigmaStep: the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
338 // pm: Decide weather Begin_Latex t > 0 End_Latex (first case) or Begin_Latex t End_Latex arbitrary (secound case)
339 // The actual work is done on the array.
341 f(x, #mu, #sigma) #Rightarrow S(t, #mu, #sigma) = (#int_{#mu}^{#mu + t #sigma} f(x, #mu, #sigma) dx + #int_{#mu}^{#mu - t #sigma} f(x, #mu, #sigma) dx) / (#int_{-#infty}^{+#infty} f(x, #mu, #sigma) dx), for t > 0
343 f(x, #mu, #sigma) #Rightarrow S(t, #mu, #sigma) = #int_{#mu}^{#mu + t #sigma} f(x, #mu, #sigma) dx / #int_{-#infty}^{+#infty} f(x, #mu, #sigma) dx
349 Float_t sigmaMax = 4;
350 gROOT->SetStyle("Plain");
351 TH1F *distribution = new TH1F("Distribution1", "Distribution f(x, #mu, #sigma)", 1000,-5,5);
353 for (Int_t i = 0; i <50000;i++) distribution->Fill(rand.Gaus(mean, sigma));
354 Float_t *ar = distribution->GetArray();
356 TCanvas* macro_example_canvas = new TCanvas("cAliBaseCalibViewer", "", 350, 350);
357 macro_example_canvas->Divide(0,3);
358 TVirtualPad *pad1 = macro_example_canvas->cd(1);
361 distribution->Draw();
362 TVirtualPad *pad2 = macro_example_canvas->cd(2);
366 TH1F *shist = AliTPCCalibViewer::SigmaCut(distribution, mean, sigma, sigmaMax);
367 shist->SetNameTitle("Cumulative","Cumulative S(t, #mu, #sigma)");
369 TVirtualPad *pad3 = macro_example_canvas->cd(3);
372 TH1F *shistPM = AliTPCCalibViewer::SigmaCut(distribution, mean, sigma, sigmaMax, -1, kTRUE);
374 return macro_example_canvas;
379 Float_t *array = histogram->GetArray();
380 Int_t nbins = histogram->GetXaxis()->GetNbins();
381 Float_t binLow = histogram->GetXaxis()->GetXmin();
382 Float_t binUp = histogram->GetXaxis()->GetXmax();
383 return AliBaseCalibViewer::SigmaCut(nbins, array, mean, sigma, nbins, binLow, binUp, sigmaMax, sigmaStep, pm);
386 //_____________________________________________________________________________
387 TH1F* AliBaseCalibViewer::SigmaCut(Int_t n, Float_t *array, Float_t mean, Float_t sigma, Int_t nbins, Float_t binLow, Float_t binUp, Float_t sigmaMax, Float_t sigmaStep, Bool_t pm){
389 // Creates a histogram Begin_Latex S(t, #mu, #sigma) End_Latex, where you can see, how much of the data are inside sigma-intervals around the mean value
390 // The data of the distribution Begin_Latex f(x, #mu, #sigma) End_Latex are given in 'array', 'n' specifies the length of the array
391 // 'mean' and 'sigma' are Begin_Latex #mu End_Latex and Begin_Latex #sigma End_Latex of the distribution in 'array', to be specified by the user
392 // 'nbins': number of bins, 'binLow': first bin, 'binUp': last bin
393 // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, Begin_Latex t #sigma End_Latex)
394 // sigmaStep: the binsize of the generated histogram
395 // Here the actual work is done.
397 if (TMath::Abs(sigma) < 1.e-10) return 0;
398 Float_t binWidth = (binUp-binLow)/(nbins - 1);
399 if (sigmaStep <= 0) sigmaStep = binWidth;
400 Int_t kbins = (Int_t)(sigmaMax * sigma / sigmaStep) + 1; // + 1 due to overflow bin in histograms
401 if (pm) kbins = 2 * (Int_t)(sigmaMax * sigma / sigmaStep) + 1;
402 Float_t kbinLow = !pm ? 0 : -sigmaMax;
403 Float_t kbinUp = sigmaMax;
404 TH1F *hist = new TH1F("sigmaCutHisto","Cumulative; Multiples of #sigma; Fraction of included data", kbins, kbinLow, kbinUp);
405 hist->SetDirectory(0);
408 // calculate normalization
409 Double_t normalization = 0;
410 for (Int_t i = 0; i <= n; i++) {
411 normalization += array[i];
414 // given units: units from given histogram
415 // sigma units: in units of sigma
416 // iDelta: integrate in interval (mean +- iDelta), given units
417 // x: ofset from mean for integration, given units
421 for (Float_t iDelta = 0; iDelta <= sigmaMax * sigma; iDelta += sigmaStep) {
423 Double_t valueP = array[GetBin(mean, nbins, binLow, binUp)];
424 Double_t valueM = array[GetBin(mean-binWidth, nbins, binLow, binUp)];
425 // add bin of mean value only once to the histogram
426 for (Float_t x = binWidth; x <= iDelta; x += binWidth) {
427 valueP += (mean + x <= binUp) ? array[GetBin(mean + x, nbins, binLow, binUp)] : 0;
428 valueM += (mean-binWidth - x >= binLow) ? array[GetBin(mean-binWidth - x, nbins, binLow, binUp)] : 0;
431 if (valueP / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail +++ \n", valueP, normalization);
432 if (valueP / normalization > 100) return hist;
433 if (valueM / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail +++ \n", valueM, normalization);
434 if (valueM / normalization > 100) return hist;
435 valueP = (valueP / normalization);
436 valueM = (valueM / normalization);
438 Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
439 hist->SetBinContent(bin, valueP);
440 bin = GetBin(-iDelta/sigma, kbins, kbinLow, kbinUp);
441 hist->SetBinContent(bin, valueM);
444 Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
445 hist->SetBinContent(bin, valueP + valueM);
448 if (!pm) hist->SetMaximum(1.2);
452 //_____________________________________________________________________________
453 TH1F* AliBaseCalibViewer::SigmaCut(Int_t n, Double_t *array, Double_t mean, Double_t sigma, Int_t nbins, Double_t *xbins, Double_t sigmaMax){
455 // SigmaCut for variable binsize
456 // NOT YET IMPLEMENTED !!!
458 printf("SigmaCut with variable binsize, Not yet implemented\n");
459 // avoid compiler warnings:
471 //_____________________________________________________________________________
472 Int_t AliBaseCalibViewer::DrawHisto1D(const Char_t* drawCommand, const Char_t* sector, const Char_t* cuts,
473 const Char_t *sigmas, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM) const
476 // Easy drawing of data, in principle the same as EasyDraw1D
477 // Difference: A line for the mean / median / LTM is drawn
478 // in 'sigmas' you can specify in which distance to the mean/median/LTM you want to see a line in sigma-units, separated by ';'
479 // example: sigmas = "2; 4; 6;" at Begin_Latex 2 #sigma End_Latex, Begin_Latex 4 #sigma End_Latex and Begin_Latex 6 #sigma End_Latex a line is drawn.
480 // "plotMean", "plotMedian" and "plotLTM": what kind of lines do you want to see?
482 Int_t oldOptStat = gStyle->GetOptStat();
483 gStyle->SetOptStat(0000000);
484 Double_t ltmFraction = 0.8;
486 TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
487 TVectorF nsigma(sigmasTokens->GetEntriesFast());
488 for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
489 TString str(((TObjString*)sigmasTokens->At(i))->GetString());
490 Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
494 TString drawStr(drawCommand);
495 Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
496 if (dangerousToDraw) {
497 Warning("DrawHisto1D", "The draw string must not contain ':' or '>>'.");
500 drawStr += " >> tempHist";
501 Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts);
502 TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
503 // FIXME is this histogram deleted automatically?
504 Double_t *values = fTree->GetV1(); // value is the array containing 'entries' numbers
506 Double_t mean = TMath::Mean(entries, values);
507 Double_t median = TMath::Median(entries, values);
508 Double_t sigma = TMath::RMS(entries, values);
509 Double_t maxY = htemp->GetMaximum();
511 TLegend * legend = new TLegend(.7,.7, .99, .99, "Statistical information");
512 //fListOfObjectsToBeDeleted->Add(legend);
516 TLine* line = new TLine(mean, 0, mean, maxY);
517 //fListOfObjectsToBeDeleted->Add(line);
518 line->SetLineColor(kRed);
519 line->SetLineWidth(2);
520 line->SetLineStyle(1);
522 legend->AddEntry(line, Form("Mean: %f", mean), "l");
524 for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
525 TLine* linePlusSigma = new TLine(mean + nsigma[i] * sigma, 0, mean + nsigma[i] * sigma, maxY);
526 //fListOfObjectsToBeDeleted->Add(linePlusSigma);
527 linePlusSigma->SetLineColor(kRed);
528 linePlusSigma->SetLineStyle(2 + i);
529 linePlusSigma->Draw();
530 TLine* lineMinusSigma = new TLine(mean - nsigma[i] * sigma, 0, mean - nsigma[i] * sigma, maxY);
531 //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
532 lineMinusSigma->SetLineColor(kRed);
533 lineMinusSigma->SetLineStyle(2 + i);
534 lineMinusSigma->Draw();
535 legend->AddEntry(lineMinusSigma, Form("%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i] * sigma)), "l");
540 TLine* line = new TLine(median, 0, median, maxY);
541 //fListOfObjectsToBeDeleted->Add(line);
542 line->SetLineColor(kBlue);
543 line->SetLineWidth(2);
544 line->SetLineStyle(1);
546 legend->AddEntry(line, Form("Median: %f", median), "l");
548 for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
549 TLine* linePlusSigma = new TLine(median + nsigma[i] * sigma, 0, median + nsigma[i]*sigma, maxY);
550 //fListOfObjectsToBeDeleted->Add(linePlusSigma);
551 linePlusSigma->SetLineColor(kBlue);
552 linePlusSigma->SetLineStyle(2 + i);
553 linePlusSigma->Draw();
554 TLine* lineMinusSigma = new TLine(median - nsigma[i] * sigma, 0, median - nsigma[i]*sigma, maxY);
555 //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
556 lineMinusSigma->SetLineColor(kBlue);
557 lineMinusSigma->SetLineStyle(2 + i);
558 lineMinusSigma->Draw();
559 legend->AddEntry(lineMinusSigma, Form("%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i] * sigma)), "l");
565 Double_t ltm = GetLTM(entries, values, <mRms, ltmFraction);
566 TLine* line = new TLine(ltm, 0, ltm, maxY);
567 //fListOfObjectsToBeDeleted->Add(line);
568 line->SetLineColor(kGreen+2);
569 line->SetLineWidth(2);
570 line->SetLineStyle(1);
572 legend->AddEntry(line, Form("LTM: %f", ltm), "l");
574 for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
575 TLine* linePlusSigma = new TLine(ltm + nsigma[i] * ltmRms, 0, ltm + nsigma[i] * ltmRms, maxY);
576 //fListOfObjectsToBeDeleted->Add(linePlusSigma);
577 linePlusSigma->SetLineColor(kGreen+2);
578 linePlusSigma->SetLineStyle(2+i);
579 linePlusSigma->Draw();
581 TLine* lineMinusSigma = new TLine(ltm - nsigma[i] * ltmRms, 0, ltm - nsigma[i] * ltmRms, maxY);
582 //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
583 lineMinusSigma->SetLineColor(kGreen+2);
584 lineMinusSigma->SetLineStyle(2+i);
585 lineMinusSigma->Draw();
586 legend->AddEntry(lineMinusSigma, Form("%i #sigma = %f", (Int_t)(nsigma[i]), (Float_t)(nsigma[i] * ltmRms)), "l");
589 if (!plotMean && !plotMedian && !plotLTM) return -1;
591 gStyle->SetOptStat(oldOptStat);
595 //_____________________________________________________________________________
596 Int_t AliBaseCalibViewer::SigmaCut(const Char_t* drawCommand, const Char_t* sector, const Char_t* cuts,
597 Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, Bool_t pm,
598 const Char_t *sigmas, Float_t sigmaStep) const {
600 // Creates a histogram, where you can see, how much of the data are inside sigma-intervals
601 // around the mean/median/LTM
602 // with drawCommand, sector and cuts you specify your input data, see EasyDraw
603 // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma)
604 // sigmaStep: the binsize of the generated histogram
605 // plotMean/plotMedian/plotLTM: specifies where to put the center
608 Double_t ltmFraction = 0.8;
610 TString drawStr(drawCommand);
611 Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
612 if (dangerousToDraw) {
613 Warning("SigmaCut", "The draw string must not contain ':' or '>>'.");
616 drawStr += " >> tempHist";
618 Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
619 TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
620 // FIXME is this histogram deleted automatically?
621 Double_t *values = fTree->GetV1(); // value is the array containing 'entries' numbers
623 Double_t mean = TMath::Mean(entries, values);
624 Double_t median = TMath::Median(entries, values);
625 Double_t sigma = TMath::RMS(entries, values);
627 TLegend * legend = new TLegend(.7,.7, .99, .99, "Cumulative");
628 //fListOfObjectsToBeDeleted->Add(legend);
629 TH1F *cutHistoMean = 0;
630 TH1F *cutHistoMedian = 0;
631 TH1F *cutHistoLTM = 0;
633 TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
634 TVectorF nsigma(sigmasTokens->GetEntriesFast());
635 for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
636 TString str(((TObjString*)sigmasTokens->At(i))->GetString());
637 Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
642 cutHistoMean = SigmaCut(htemp, mean, sigma, sigmaMax, sigmaStep, pm);
644 //fListOfObjectsToBeDeleted->Add(cutHistoMean);
645 cutHistoMean->SetLineColor(kRed);
646 legend->AddEntry(cutHistoMean, "Mean", "l");
647 cutHistoMean->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
648 cutHistoMean->Draw();
649 DrawLines(cutHistoMean, nsigma, legend, kRed, pm);
650 } // if (cutHistoMean)
654 cutHistoMedian = SigmaCut(htemp, median, sigma, sigmaMax, sigmaStep, pm);
655 if (cutHistoMedian) {
656 //fListOfObjectsToBeDeleted->Add(cutHistoMedian);
657 cutHistoMedian->SetLineColor(kBlue);
658 legend->AddEntry(cutHistoMedian, "Median", "l");
659 cutHistoMedian->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
660 if (plotMean && cutHistoMean) cutHistoMedian->Draw("same");
661 else cutHistoMedian->Draw();
662 DrawLines(cutHistoMedian, nsigma, legend, kBlue, pm);
663 } // if (cutHistoMedian)
667 Double_t ltm = GetLTM(entries, values, <mRms, ltmFraction);
668 cutHistoLTM = SigmaCut(htemp, ltm, ltmRms, sigmaMax, sigmaStep, pm);
670 //fListOfObjectsToBeDeleted->Add(cutHistoLTM);
671 cutHistoLTM->SetLineColor(kGreen+2);
672 legend->AddEntry(cutHistoLTM, "LTM", "l");
673 cutHistoLTM->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
674 if ((plotMean && cutHistoMean) || (plotMedian && cutHistoMedian)) cutHistoLTM->Draw("same");
675 else cutHistoLTM->Draw();
676 DrawLines(cutHistoLTM, nsigma, legend, kGreen+2, pm);
679 if (!plotMean && !plotMedian && !plotLTM) return -1;
684 //_____________________________________________________________________________
685 Int_t AliBaseCalibViewer::Integrate(const Char_t* drawCommand, const Char_t* sector, const Char_t* cuts,
686 Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM,
687 const Char_t *sigmas, Float_t sigmaStep) const {
689 // Creates an integrated histogram Begin_Latex S(t, #mu, #sigma) End_Latex, out of the input distribution distribution Begin_Latex f(x, #mu, #sigma) End_Latex, given in "histogram"
690 // "mean" and "sigma" are Begin_Latex #mu End_Latex and Begin_Latex #sigma End_Latex of the distribution in "histogram", to be specified by the user
691 // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate
692 // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
693 // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
694 // The actual work is done on the array.
696 f(x, #mu, #sigma) #Rightarrow S(t, #mu, #sigma) = #int_{-#infty}^{#mu + t #sigma} f(x, #mu, #sigma) dx / #int_{-#infty}^{+#infty} f(x, #mu, #sigma) dx
700 Double_t ltmFraction = 0.8;
701 // avoid compiler warnings:
703 sigmaStep = sigmaStep;
705 TString drawStr(drawCommand);
706 Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
707 if (dangerousToDraw) {
708 Warning("Integrate", "The draw string must not contain ':' or '>>'.");
711 drawStr += " >> tempHist";
713 Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
714 TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
715 TGraph *integralGraphMean = 0;
716 TGraph *integralGraphMedian = 0;
717 TGraph *integralGraphLTM = 0;
718 Double_t *values = fTree->GetV1(); // value is the array containing 'entries' numbers
719 Int_t *index = new Int_t[entries];
720 Float_t *xarray = new Float_t[entries];
721 Float_t *yarray = new Float_t[entries];
722 TMath::Sort(entries, values, index, kFALSE);
724 Double_t mean = TMath::Mean(entries, values);
725 Double_t median = TMath::Median(entries, values);
726 Double_t sigma = TMath::RMS(entries, values);
728 // parse sigmas string
729 TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
730 TVectorF nsigma(sigmasTokens->GetEntriesFast());
731 for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
732 TString str(((TObjString*)sigmasTokens->At(i))->GetString());
733 Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
737 TLegend * legend = new TLegend(.7,.7, .99, .99, "Integrated histogram");
738 //fListOfObjectsToBeDeleted->Add(legend);
741 for (Int_t i = 0; i < entries; i++) {
742 xarray[i] = (values[index[i]] - mean) / sigma;
743 yarray[i] = float(i) / float(entries);
745 integralGraphMean = new TGraph(entries, xarray, yarray);
746 if (integralGraphMean) {
747 //fListOfObjectsToBeDeleted->Add(integralGraphMean);
748 integralGraphMean->SetLineColor(kRed);
749 legend->AddEntry(integralGraphMean, "Mean", "l");
750 integralGraphMean->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
751 integralGraphMean->Draw("alu");
752 DrawLines(integralGraphMean, nsigma, legend, kRed, kTRUE);
756 for (Int_t i = 0; i < entries; i++) {
757 xarray[i] = (values[index[i]] - median) / sigma;
758 yarray[i] = float(i) / float(entries);
760 integralGraphMedian = new TGraph(entries, xarray, yarray);
761 if (integralGraphMedian) {
762 //fListOfObjectsToBeDeleted->Add(integralGraphMedian);
763 integralGraphMedian->SetLineColor(kBlue);
764 legend->AddEntry(integralGraphMedian, "Median", "l");
765 integralGraphMedian->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
766 if (plotMean && integralGraphMean) integralGraphMedian->Draw("samelu");
767 else integralGraphMedian->Draw("alu");
768 DrawLines(integralGraphMedian, nsigma, legend, kBlue, kTRUE);
773 Double_t ltm = GetLTM(entries, values, <mRms, ltmFraction);
774 for (Int_t i = 0; i < entries; i++) {
775 xarray[i] = (values[index[i]] - ltm) / ltmRms;
776 yarray[i] = float(i) / float(entries);
778 integralGraphLTM = new TGraph(entries, xarray, yarray);
779 if (integralGraphLTM) {
780 //fListOfObjectsToBeDeleted->Add(integralGraphLTM);
781 integralGraphLTM->SetLineColor(kGreen+2);
782 legend->AddEntry(integralGraphLTM, "LTM", "l");
783 integralGraphLTM->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
784 if ((plotMean && integralGraphMean) || (plotMedian && integralGraphMedian)) integralGraphLTM->Draw("samelu");
785 else integralGraphLTM->Draw("alu");
786 DrawLines(integralGraphLTM, nsigma, legend, kGreen+2, kTRUE);
792 if (!plotMean && !plotMedian && !plotLTM) return -1;
797 //_____________________________________________________________________________
798 TH1F* AliBaseCalibViewer::Integrate(TH1F *histogram, Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep){
800 // Creates an integrated histogram Begin_Latex S(t, #mu, #sigma) End_Latex, out of the input distribution distribution Begin_Latex f(x, #mu, #sigma) End_Latex, given in "histogram"
801 // "mean" and "sigma" are Begin_Latex #mu End_Latex and Begin_Latex #sigma End_Latex of the distribution in "histogram", to be specified by the user
802 // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate
803 // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
804 // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
805 // The actual work is done on the array.
807 f(x, #mu, #sigma) #Rightarrow S(t, #mu, #sigma) = #int_{-#infty}^{#mu + t #sigma} f(x, #mu, #sigma) dx / #int_{-#infty}^{+#infty} f(x, #mu, #sigma) dx
813 Float_t sigmaMax = 4;
814 gROOT->SetStyle("Plain");
815 TH1F *distribution = new TH1F("Distribution2", "Distribution f(x, #mu, #sigma)", 1000,-5,5);
817 for (Int_t i = 0; i <50000;i++) distribution->Fill(rand.Gaus(mean, sigma));
818 Float_t *ar = distribution->GetArray();
820 TCanvas* macro_example_canvas = new TCanvas("macro_example_canvas_Integrate", "", 350, 350);
821 macro_example_canvas->Divide(0,2);
822 TVirtualPad *pad1 = macro_example_canvas->cd(1);
825 distribution->Draw();
826 TVirtualPad *pad2 = macro_example_canvas->cd(2);
829 TH1F *shist = AliTPCCalibViewer::Integrate(distribution, mean, sigma, sigmaMax);
830 shist->SetNameTitle("Cumulative","Cumulative S(t, #mu, #sigma)");
833 return macro_example_canvas_Integrate;
839 Float_t *array = histogram->GetArray();
840 Int_t nbins = histogram->GetXaxis()->GetNbins();
841 Float_t binLow = histogram->GetXaxis()->GetXmin();
842 Float_t binUp = histogram->GetXaxis()->GetXmax();
843 return Integrate(nbins, array, nbins, binLow, binUp, mean, sigma, sigmaMax, sigmaStep);
846 //_____________________________________________________________________________
847 TH1F* AliBaseCalibViewer::Integrate(Int_t n, Float_t *array, Int_t nbins, Float_t binLow, Float_t binUp,
848 Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep){
849 // Creates an integrated histogram Begin_Latex S(t, #mu, #sigma) End_Latex, out of the input distribution distribution Begin_Latex f(x, #mu, #sigma) End_Latex, given in "histogram"
850 // "mean" and "sigma" are Begin_Latex #mu End_Latex and Begin_Latex #sigma End_Latex of the distribution in "histogram", to be specified by the user
851 // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate
852 // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
853 // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
854 // Here the actual work is done.
856 Bool_t givenUnits = kTRUE;
857 if (TMath::Abs(sigma) < 1.e-10 && TMath::Abs(sigmaMax) < 1.e-10) givenUnits = kFALSE;
860 sigmaMax = (binUp - binLow) / 2.;
863 Float_t binWidth = (binUp-binLow)/(nbins - 1);
864 if (sigmaStep <= 0) sigmaStep = binWidth;
865 Int_t kbins = (Int_t)(sigmaMax * sigma / sigmaStep) + 1; // + 1 due to overflow bin in histograms
866 Float_t kbinLow = givenUnits ? binLow : -sigmaMax;
867 Float_t kbinUp = givenUnits ? binUp : sigmaMax;
869 if (givenUnits) hist = new TH1F("integratedHisto","Integrated Histogram; Given x; Fraction of included data", kbins, kbinLow, kbinUp);
870 if (!givenUnits) hist = new TH1F("integratedHisto","Integrated Histogram; Multiples of #sigma; Fraction of included data", kbins, kbinLow, kbinUp);
871 hist->SetDirectory(0);
874 // calculate normalization
875 // printf("calculating normalization, integrating from bin 1 to %i \n", n);
876 Double_t normalization = 0;
877 for (Int_t i = 1; i <= n; i++) {
878 normalization += array[i];
880 // printf("normalization: %f \n", normalization);
882 // given units: units from given histogram
883 // sigma units: in units of sigma
884 // iDelta: integrate in interval (mean +- iDelta), given units
885 // x: ofset from mean for integration, given units
889 for (Float_t iDelta = mean - sigmaMax * sigma; iDelta <= mean + sigmaMax * sigma; iDelta += sigmaStep) {
892 for (Float_t x = mean - sigmaMax * sigma; x <= iDelta; x += binWidth) {
893 value += (x <= binUp && x >= binLow) ? array[GetBin(x, nbins, binLow, binUp)] : 0;
895 if (value / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail +++ \n", value, normalization);
896 if (value / normalization > 100) return hist;
897 Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
898 // printf("first integration bin: %i, last integration bin: %i \n", GetBin(mean - sigmaMax * sigma, nbins, binLow, binUp), GetBin(iDelta, nbins, binLow, binUp));
899 // printf("value: %f, normalization: %f, normalized value: %f, iDelta: %f, Bin: %i \n", value, normalization, value/normalization, iDelta, bin);
900 value = (value / normalization);
901 hist->SetBinContent(bin, value);
906 //_____________________________________________________________________________
907 void AliBaseCalibViewer::DrawLines(TH1F *histogram, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const {
909 // Private function for SigmaCut(...) and Integrate(...)
910 // Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
913 // start to draw the lines, loop over requested sigmas
914 for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
916 Int_t bin = histogram->GetXaxis()->FindBin(nsigma[i]);
917 TLine* lineUp = new TLine(nsigma[i], 0, nsigma[i], histogram->GetBinContent(bin));
918 //fListOfObjectsToBeDeleted->Add(lineUp);
919 lineUp->SetLineColor(color);
920 lineUp->SetLineStyle(2 + i);
922 TLine* lineLeft = new TLine(nsigma[i], histogram->GetBinContent(bin), 0, histogram->GetBinContent(bin));
923 //fListOfObjectsToBeDeleted->Add(lineLeft);
924 lineLeft->SetLineColor(color);
925 lineLeft->SetLineStyle(2 + i);
927 legend->AddEntry(lineLeft, Form("Fraction(%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin)), "l");
930 Int_t bin = histogram->GetXaxis()->FindBin(nsigma[i]);
931 TLine* lineUp1 = new TLine(nsigma[i], 0, nsigma[i], histogram->GetBinContent(bin));
932 //fListOfObjectsToBeDeleted->Add(lineUp1);
933 lineUp1->SetLineColor(color);
934 lineUp1->SetLineStyle(2 + i);
936 TLine* lineLeft1 = new TLine(nsigma[i], histogram->GetBinContent(bin), histogram->GetBinLowEdge(0)+histogram->GetBinWidth(0), histogram->GetBinContent(bin));
937 //fListOfObjectsToBeDeleted->Add(lineLeft1);
938 lineLeft1->SetLineColor(color);
939 lineLeft1->SetLineStyle(2 + i);
941 legend->AddEntry(lineLeft1, Form("Fraction(+%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin)), "l");
942 bin = histogram->GetXaxis()->FindBin(-nsigma[i]);
943 TLine* lineUp2 = new TLine(-nsigma[i], 0, -nsigma[i], histogram->GetBinContent(bin));
944 //fListOfObjectsToBeDeleted->Add(lineUp2);
945 lineUp2->SetLineColor(color);
946 lineUp2->SetLineStyle(2 + i);
948 TLine* lineLeft2 = new TLine(-nsigma[i], histogram->GetBinContent(bin), histogram->GetBinLowEdge(0)+histogram->GetBinWidth(0), histogram->GetBinContent(bin));
949 //fListOfObjectsToBeDeleted->Add(lineLeft2);
950 lineLeft2->SetLineColor(color);
951 lineLeft2->SetLineStyle(2 + i);
953 legend->AddEntry(lineLeft2, Form("Fraction(-%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin)), "l");
955 } // for (Int_t i = 0; i < nsigma.GetNoElements(); i++)
958 //_____________________________________________________________________________
959 void AliBaseCalibViewer::DrawLines(TGraph *graph, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const {
961 // Private function for SigmaCut(...) and Integrate(...)
962 // Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
965 // start to draw the lines, loop over requested sigmas
966 for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
968 TLine* lineUp = new TLine(nsigma[i], 0, nsigma[i], graph->Eval(nsigma[i]));
969 //fListOfObjectsToBeDeleted->Add(lineUp);
970 lineUp->SetLineColor(color);
971 lineUp->SetLineStyle(2 + i);
973 TLine* lineLeft = new TLine(nsigma[i], graph->Eval(nsigma[i]), 0, graph->Eval(nsigma[i]));
974 //fListOfObjectsToBeDeleted->Add(lineLeft);
975 lineLeft->SetLineColor(color);
976 lineLeft->SetLineStyle(2 + i);
978 legend->AddEntry(lineLeft, Form("Fraction(%f #sigma) = %f",nsigma[i], graph->Eval(nsigma[i])), "l");
981 TLine* lineUp1 = new TLine(nsigma[i], 0, nsigma[i], graph->Eval(nsigma[i]));
982 //fListOfObjectsToBeDeleted->Add(lineUp1);
983 lineUp1->SetLineColor(color);
984 lineUp1->SetLineStyle(2 + i);
986 TLine* lineLeft1 = new TLine(nsigma[i], graph->Eval(nsigma[i]), graph->GetHistogram()->GetXaxis()->GetBinLowEdge(0), graph->Eval(nsigma[i]));
987 //fListOfObjectsToBeDeleted->Add(lineLeft1);
988 lineLeft1->SetLineColor(color);
989 lineLeft1->SetLineStyle(2 + i);
991 legend->AddEntry(lineLeft1, Form("Fraction(+%f #sigma) = %f",nsigma[i], graph->Eval(nsigma[i])), "l");
992 TLine* lineUp2 = new TLine(-nsigma[i], 0, -nsigma[i], graph->Eval(-nsigma[i]));
993 //fListOfObjectsToBeDeleted->Add(lineUp2);
994 lineUp2->SetLineColor(color);
995 lineUp2->SetLineStyle(2 + i);
997 TLine* lineLeft2 = new TLine(-nsigma[i], graph->Eval(-nsigma[i]), graph->GetHistogram()->GetXaxis()->GetBinLowEdge(0), graph->Eval(-nsigma[i]));
998 //fListOfObjectsToBeDeleted->Add(lineLeft2);
999 lineLeft2->SetLineColor(color);
1000 lineLeft2->SetLineStyle(2 + i);
1002 legend->AddEntry(lineLeft2, Form("Fraction(-%f #sigma) = %f",nsigma[i], graph->Eval(-nsigma[i])), "l");
1004 } // for (Int_t i = 0; i < nsigma.GetNoElements(); i++)