1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
17 ///////////////////////////////////////////////////////////////////////////////
19 // Class for viewing/visualizing TPC calibration data //
20 // base on TTree functionality for visualization //
22 // Create a list of AliTPCCalPads, arrange them in an TObjArray. //
23 // Pass this TObjArray to MakeTree and create the calibration Tree //
24 // While craating this tree some statistical information are calculated //
25 // Open the viewer with this Tree: AliTPCCalibViewer v("CalibTree.root") //
27 // EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0") //
29 // If you like to click, we recommand you the //
30 // AliTPCCalibViewerGUI //
32 // THE DOCUMENTATION IS STILL NOT COMPLETED !!!! //
34 ///////////////////////////////////////////////////////////////////////////////
50 #include <THashTable.h>
51 #include <TObjString.h>
52 #include "TTreeStream.h"
56 #include "TDirectory.h"
57 #include "AliTPCCalibPulser.h"
58 #include "AliTPCCalibPedestal.h"
59 #include "AliTPCCalibCE.h"
60 #include "TFriendElement.h"
61 // #include "TObjArray.h"
62 // #include "TObjString.h"
63 // #include "TString.h"
64 // #include "AliTPCCalPad.h"
70 #include "AliTPCCalibViewer.h"
72 ClassImp(AliTPCCalibViewer)
75 AliTPCCalibViewer::AliTPCCalibViewer()
79 fListOfObjectsToBeDeleted(0),
80 fTreeMustBeDeleted(0),
85 // Default constructor
90 //_____________________________________________________________________________
91 AliTPCCalibViewer::AliTPCCalibViewer(const AliTPCCalibViewer &c)
95 fListOfObjectsToBeDeleted(0),
96 fTreeMustBeDeleted(0),
101 // dummy AliTPCCalibViewer copy constructor
102 // not yet working!!!
105 fTreeMustBeDeleted = c.fTreeMustBeDeleted;
106 //fFile = new TFile(*(c.fFile));
107 fListOfObjectsToBeDeleted = c.fListOfObjectsToBeDeleted;
108 fAbbreviation = c.fAbbreviation;
109 fAppendString = c.fAppendString;
112 //_____________________________________________________________________________
113 AliTPCCalibViewer::AliTPCCalibViewer(TTree* tree)
117 fListOfObjectsToBeDeleted(0),
118 fTreeMustBeDeleted(0),
123 // Constructor that initializes the calibration viewer
126 fTreeMustBeDeleted = kFALSE;
127 fListOfObjectsToBeDeleted = new TObjArray();
129 fAppendString = ".fElements";
132 //_____________________________________________________________________________
133 AliTPCCalibViewer::AliTPCCalibViewer(char* fileName, char* treeName)
137 fListOfObjectsToBeDeleted(0),
138 fTreeMustBeDeleted(0),
144 // Constructor to initialize the calibration viewer
145 // the file 'fileName' contains the tree 'treeName'
147 fFile = new TFile(fileName, "read");
148 fTree = (TTree*) fFile->Get(treeName);
149 fTreeMustBeDeleted = kTRUE;
150 fListOfObjectsToBeDeleted = new TObjArray();
152 fAppendString = ".fElements";
155 //____________________________________________________________________________
156 AliTPCCalibViewer & AliTPCCalibViewer::operator =(const AliTPCCalibViewer & param)
159 // assignment operator - dummy
160 // not yet working!!!
163 fTreeMustBeDeleted = param.fTreeMustBeDeleted;
164 //fFile = new TFile(*(param.fFile));
165 fListOfObjectsToBeDeleted = param.fListOfObjectsToBeDeleted;
166 fAbbreviation = param.fAbbreviation;
167 fAppendString = param.fAppendString;
171 //_____________________________________________________________________________
172 AliTPCCalibViewer::~AliTPCCalibViewer()
175 // AliTPCCalibViewer destructor
176 // all objects will be deleted, the file will be closed, the pictures will disappear
178 if (fTree && fTreeMustBeDeleted) {
179 fTree->SetCacheSize(0);
188 for (Int_t i = fListOfObjectsToBeDeleted->GetEntriesFast()-1; i >= 0; i--) {
189 //cout << "Index " << i << " trying to delete the following object: " << fListOfObjectsToBeDeleted->At(i)->GetName() << "..."<< endl;
190 delete fListOfObjectsToBeDeleted->At(i);
192 delete fListOfObjectsToBeDeleted;
195 //_____________________________________________________________________________
196 void AliTPCCalibViewer::Delete(Option_t* option) {
198 // Should be called from AliTPCCalibViewerGUI class only.
199 // If you use Delete() do not call the destructor.
200 // All objects (except those contained in fListOfObjectsToBeDeleted) will be deleted, the file will be closed.
203 option = option; // to avoid warnings on compiling
204 if (fTree && fTreeMustBeDeleted) {
205 fTree->SetCacheSize(0);
210 delete fListOfObjectsToBeDeleted;
214 const char* AliTPCCalibViewer::AddAbbreviations(char* c, Bool_t printDrawCommand){
215 // Replace all "<variable>" with "<variable><fAbbreviation>" (Adds forgotten "~")
216 // but take care on the statistical information, like "CEQmean_Mean"
217 // and also take care on correct given variables, like "CEQmean~"
219 // For each variable out of "listOfVariables":
220 // - 'Save' correct items:
221 // - form <replaceString>, take <variable>'s first char, add <removeString>, add rest of <variable>, e.g. "C!#EQmean" (<removeString> = "!#")
222 // - For each statistical information in "listOfNormalizationVariables":
223 // - ReplaceAll <variable><statistical_Information> with <replaceString><statistical_Information>
224 // - ReplaceAll <variable><abbreviation> with <replaceString><abbreviation>, e.g. "CEQmean~" -> "C!#EQmean~"
225 // - ReplaceAll <variable><appendStr> with <replaceString><appendStr>, e.g. "CEQmean.fElements" -> "C!#EQmean.fElements"
227 // - Do actual replacing:
228 // - ReplaceAll <variable> with <variable><fAbbreviation>, e.g. "CEQmean" -> "CEQmean~"
231 // - For each statistical information in "listOfNormalizationVariables":
232 // - ReplaceAll <replaceString><statistical_Information> with <variable><statistical_Information>
233 // - ReplaceAll <replaceString><abbreviation> with <variable><abbreviation>, e.g. "C!#EQmean~" -> "CEQmean~"
234 // - ReplaceAll <replaceString><appendStr> with <variable><appendStr>, e.g. "C!#EQmean.fElements" -> "CEQmean.fElements"
236 // Now all the missing "~" should be added.
239 TString removeString = "!#"; // very unpropable combination of chars
240 TString replaceString = "";
241 TString searchString = "";
242 TString normString = "";
243 TObjArray *listOfVariables = GetListOfVariables();
244 listOfVariables->Add(new TObjString("channel"));
245 listOfVariables->Add(new TObjString("gx"));
246 listOfVariables->Add(new TObjString("gy"));
247 listOfVariables->Add(new TObjString("lx"));
248 listOfVariables->Add(new TObjString("ly"));
249 listOfVariables->Add(new TObjString("pad"));
250 listOfVariables->Add(new TObjString("row"));
251 listOfVariables->Add(new TObjString("rpad"));
252 listOfVariables->Add(new TObjString("sector"));
253 TObjArray *listOfNormalizationVariables = GetListOfNormalizationVariables();
254 Int_t nVariables = listOfVariables->GetEntriesFast();
255 Int_t nNorm = listOfNormalizationVariables->GetEntriesFast();
257 Int_t *varLengths = new Int_t[nVariables];
258 for (Int_t i = 0; i < nVariables; i++) {
259 varLengths[i] = ((TObjString*)listOfVariables->At(i))->String().Length();
261 Int_t *normLengths = new Int_t[nNorm];
262 for (Int_t i = 0; i < nNorm; i++) {
263 normLengths[i] = ((TObjString*)listOfNormalizationVariables->At(i))->String().Length();
264 // printf("normLengths[%i] (%s) = %i \n", i,((TObjString*)listOfNormalizationVariables->At(i))->String().Data(), normLengths[i]);
266 Int_t *varSort = new Int_t[nVariables];
267 TMath::Sort(nVariables, varLengths, varSort, kTRUE);
268 Int_t *normSort = new Int_t[nNorm];
269 TMath::Sort(nNorm, normLengths, normSort, kTRUE);
270 // for (Int_t i = 0; i<nNorm; i++) printf("normLengths: %i\n", normLengths[normSort[i]]);
271 // for (Int_t i = 0; i<nVariables; i++) printf("varLengths: %i\n", varLengths[varSort[i]]);
273 for (Int_t ivar = 0; ivar < nVariables; ivar++) {
274 // ***** save correct tokens *****
275 // first get the next variable:
276 searchString = ((TObjString*)listOfVariables->At(varSort[ivar]))->String();
277 // printf("searchString: %s ++++++++++++++\n", searchString.Data());
278 // form replaceString:
280 for (Int_t i = 0; i < searchString.Length(); i++) {
281 replaceString.Append(searchString[i]);
282 if (i == 0) replaceString.Append(removeString);
284 // go through normalization:
285 // printf("go through normalization\n");
286 for (Int_t inorm = 0; inorm < nNorm; inorm++) {
287 // printf(" inorm=%i, nNorm=%i, normSort[inorm]=%i \n", inorm, nNorm, normSort[inorm]);
288 normString = ((TObjString*)listOfNormalizationVariables->At(normSort[inorm]))->String();
289 // printf(" walking in normalization, i=%i, normString=%s \n", inorm, normString.Data());
290 str.ReplaceAll(searchString + normString, replaceString + normString);
291 // like: str.ReplaceAll("CEQmean_Mean", "C!EQmean_Mean");
293 str.ReplaceAll(searchString + fAbbreviation, replaceString + fAbbreviation);
294 // like: str.ReplaceAll("CEQmean~", "C!EQmean~");
295 str.ReplaceAll(searchString + fAppendString, replaceString + fAppendString);
296 // like: str.ReplaceAll("CEQmean.fElements", "C!EQmean.fElements");
298 // ***** add missing extensions *****
299 str.ReplaceAll(searchString, replaceString + fAbbreviation);
300 // like: str.ReplaceAll("CEQmean", "C!EQmean~");
303 // ***** undo saving *****
304 str.ReplaceAll(removeString, "");
306 if (printDrawCommand) std::cout << "The string looks now like: " << str.Data() << std::endl;
315 //_____________________________________________________________________________
316 Int_t AliTPCCalibViewer::EasyDraw(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
318 // easy drawing of data, use '~' for abbreviation of '.fElements'
319 // example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
320 // sector: sector-number - only the specified sector will be drwawn
321 // 'A'/'C' or 'a'/'c' - side A/C will be drawn
322 // 'ALL' - whole TPC will be drawn, projected on one side
323 // cuts: specifies cuts
324 // drawOptions: draw options like 'same'
325 // writeDrawCommand: write the command, that is passed to TTree::Draw
328 TString drawStr(drawCommand);
329 TString sectorStr(sector);
332 //TString drawOptionsStr("profcolz ");
333 Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
334 if (dangerousToDraw) {
335 Warning("EasyDraw", "The draw string must not contain ':' or '>>'.");
338 TString drawOptionsStr("");
340 Int_t rndNumber = rnd.Integer(10000);
342 if (drawOptions && strcmp(drawOptions, "") != 0)
343 drawOptionsStr += drawOptions;
345 drawOptionsStr += "profcolz";
347 if (sectorStr == "A") {
348 drawStr += Form(":gy%s:gx%s>>prof", fAppendString.Data(), fAppendString.Data());
349 drawStr += rndNumber;
350 drawStr += "(330,-250,250,330,-250,250)";
351 cutStr += "(sector/18)%2==0 ";
353 else if (sectorStr == "C") {
354 drawStr += Form(":gy%s:gx%s>>prof", fAppendString.Data(), fAppendString.Data());
355 drawStr += rndNumber;
356 drawStr += "(330,-250,250,330,-250,250)";
357 cutStr += "(sector/18)%2==1 ";
359 else if (sectorStr == "ALL") {
360 drawStr += Form(":gy%s:gx%s>>prof", fAppendString.Data(), fAppendString.Data());
361 drawStr += rndNumber;
362 drawStr += "(330,-250,250,330,-250,250)";
364 else if (sectorStr.IsDigit()) {
365 Int_t isec = sectorStr.Atoi();
366 drawStr += Form(":rpad%s:row%s>>prof", fAppendString.Data(), fAppendString.Data());
367 drawStr += rndNumber;
368 if (isec < 36 && isec >= 0)
369 drawStr += "(63,0,63,108,-54,54)";
370 else if (isec < 72 && isec >= 36)
371 drawStr += "(96,0,96,140,-70,70)";
373 Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
376 cutStr += "(sector==";
381 if (cuts && cuts[0] != 0) {
382 if (cutStr.Length() != 0) cutStr += "&& ";
387 drawStr.ReplaceAll(fAbbreviation, fAppendString);
388 cutStr.ReplaceAll(fAbbreviation, fAppendString);
389 if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" << cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
390 Int_t returnValue = fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
391 TString profName("prof");
392 profName += rndNumber;
393 TObject *obj = gDirectory->Get(profName.Data());
394 if (obj && obj->InheritsFrom("TH1")) FormatHistoLabels((TH1*)obj);
399 Int_t AliTPCCalibViewer::EasyDraw(const char* drawCommand, Int_t sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
401 // easy drawing of data, use '~' for abbreviation of '.fElements'
402 // example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
403 // sector: sector-number - only the specified sector will be drwawn
404 // cuts: specifies cuts
405 // drawOptions: draw options like 'same'
406 // writeDrawCommand: write the command, that is passed to TTree::Draw
408 if (sector >= 0 && sector < 72) {
410 sprintf(sectorChr, "%i", sector);
411 return EasyDraw(drawCommand, sectorChr, cuts, drawOptions, writeDrawCommand);
413 Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
418 //_____________________________________________________________________________
419 Int_t AliTPCCalibViewer::EasyDraw1D(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
421 // easy drawing of data, use '~' for abbreviation of '.fElements'
422 // example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
423 // sector: sector-number - the specified sector will be drwawn
424 // 'A'/'C' or 'a'/'c' - side A/C will be drawn
425 // 'ALL' - whole TPC will be drawn, projected on one side
426 // cuts: specifies cuts
427 // drawOptions: draw options like 'same'
428 // writeDrawCommand: write the command, that is passed to TTree::Draw
431 TString drawStr(drawCommand);
432 TString sectorStr(sector);
433 TString drawOptionsStr(drawOptions);
437 if (sectorStr == "A")
438 cutStr += "(sector/18)%2==0 ";
439 else if (sectorStr == "C")
440 cutStr += "(sector/18)%2==1 ";
441 else if (sectorStr.IsDigit()) {
442 Int_t isec = sectorStr.Atoi();
443 if (isec < 0 || isec > 71) {
444 Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
447 cutStr += "(sector==";
452 if (cuts && cuts[0] != 0) {
453 if (cutStr.Length() != 0) cutStr += "&& ";
459 drawStr.ReplaceAll(fAbbreviation, fAppendString);
460 cutStr.ReplaceAll(fAbbreviation, fAppendString);
461 if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" << cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
462 Int_t returnValue = fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
463 if (returnValue == -1) return -1;
465 TObject *obj = (gPad) ? gPad->GetPrimitive("htemp") : 0;
466 if (!obj) obj = (TH1F*)gDirectory->Get("htemp");
467 if (!obj) obj = gPad->GetPrimitive("tempHist");
468 if (!obj) obj = (TH1F*)gDirectory->Get("tempHist");
469 if (!obj) obj = gPad->GetPrimitive("Graph");
470 if (!obj) obj = (TH1F*)gDirectory->Get("Graph");
471 if (obj && obj->InheritsFrom("TH1")) FormatHistoLabels((TH1*)obj);
476 Int_t AliTPCCalibViewer::EasyDraw1D(const char* drawCommand, Int_t sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
478 // easy drawing of data, use '~' for abbreviation of '.fElements'
479 // example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
480 // sector: sector-number - the specified sector will be drwawn
481 // cuts: specifies cuts
482 // drawOptions: draw options like 'same'
483 // writeDrawCommand: write the command, that is passed to TTree::Draw
486 if (sector >= 0 && sector < 72) {
488 sprintf(sectorChr, "%i", sector);
489 return EasyDraw1D(drawCommand, sectorChr, cuts, drawOptions, writeDrawCommand);
491 Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
496 void AliTPCCalibViewer::FormatHistoLabels(TH1 *histo) const {
498 // formats title and axis labels of histo
499 // removes '.fElements'
502 TString replaceString(fAppendString.Data());
503 TString *str = new TString(histo->GetTitle());
504 str->ReplaceAll(replaceString, "");
505 histo->SetTitle(str->Data());
507 if (histo->GetXaxis()) {
508 str = new TString(histo->GetXaxis()->GetTitle());
509 str->ReplaceAll(replaceString, "");
510 histo->GetXaxis()->SetTitle(str->Data());
513 if (histo->GetYaxis()) {
514 str = new TString(histo->GetYaxis()->GetTitle());
515 str->ReplaceAll(replaceString, "");
516 histo->GetYaxis()->SetTitle(str->Data());
519 if (histo->GetZaxis()) {
520 str = new TString(histo->GetZaxis()->GetTitle());
521 str->ReplaceAll(replaceString, "");
522 histo->GetZaxis()->SetTitle(str->Data());
528 Int_t AliTPCCalibViewer::DrawHisto1D(const char* drawCommand, Int_t sector, const char* cuts, const char *sigmas, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM) const {
530 // Easy drawing of data, in principle the same as EasyDraw1D
531 // Difference: A line for the mean / median / LTM is drawn
532 // in 'sigmas' you can specify in which distance to the mean/median/LTM you want to see a line in sigma-units, separated by ';'
533 // 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.
534 // "plotMean", "plotMedian" and "plotLTM": what kind of lines do you want to see?
536 if (sector >= 0 && sector < 72) {
538 sprintf(sectorChr, "%i", sector);
539 return DrawHisto1D(drawCommand, sectorChr, cuts, sigmas, plotMean, plotMedian, plotLTM);
541 Error("DrawHisto1D","The TPC contains only sectors between 0 and 71.");
546 Int_t AliTPCCalibViewer::DrawHisto1D(const char* drawCommand, const char* sector, const char* cuts, const char *sigmas, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM) const {
548 // Easy drawing of data, in principle the same as EasyDraw1D
549 // Difference: A line for the mean / median / LTM is drawn
550 // in 'sigmas' you can specify in which distance to the mean/median/LTM you want to see a line in sigma-units, separated by ';'
551 // 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.
552 // "plotMean", "plotMedian" and "plotLTM": what kind of lines do you want to see?
554 Int_t oldOptStat = gStyle->GetOptStat();
555 gStyle->SetOptStat(0000000);
556 Double_t ltmFraction = 0.8;
558 TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
559 TVectorF nsigma(sigmasTokens->GetEntriesFast());
560 for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
561 TString str(((TObjString*)sigmasTokens->At(i))->GetString());
562 Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
566 TString drawStr(drawCommand);
567 Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
568 if (dangerousToDraw) {
569 Warning("DrawHisto1D", "The draw string must not contain ':' or '>>'.");
572 drawStr += " >> tempHist";
573 Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts);
574 TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
575 // FIXME is this histogram deleted automatically?
576 Double_t *values = fTree->GetV1(); // value is the array containing 'entries' numbers
578 Double_t mean = TMath::Mean(entries, values);
579 Double_t median = TMath::Median(entries, values);
580 Double_t sigma = TMath::RMS(entries, values);
581 Double_t maxY = htemp->GetMaximum();
584 TLegend * legend = new TLegend(.7,.7, .99, .99, "Statistical information");
585 // sprintf(c, "%s, sector: %i", type, sector);
586 //fListOfObjectsToBeDeleted->Add(legend);
590 TLine* line = new TLine(mean, 0, mean, maxY);
591 //fListOfObjectsToBeDeleted->Add(line);
592 line->SetLineColor(kRed);
593 line->SetLineWidth(2);
594 line->SetLineStyle(1);
596 sprintf(c, "Mean: %f", mean);
597 legend->AddEntry(line, c, "l");
599 for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
600 TLine* linePlusSigma = new TLine(mean + nsigma[i] * sigma, 0, mean + nsigma[i] * sigma, maxY);
601 //fListOfObjectsToBeDeleted->Add(linePlusSigma);
602 linePlusSigma->SetLineColor(kRed);
603 linePlusSigma->SetLineStyle(2 + i);
604 linePlusSigma->Draw();
605 TLine* lineMinusSigma = new TLine(mean - nsigma[i] * sigma, 0, mean - nsigma[i] * sigma, maxY);
606 //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
607 lineMinusSigma->SetLineColor(kRed);
608 lineMinusSigma->SetLineStyle(2 + i);
609 lineMinusSigma->Draw();
610 sprintf(c, "%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i] * sigma));
611 legend->AddEntry(lineMinusSigma, c, "l");
616 TLine* line = new TLine(median, 0, median, maxY);
617 //fListOfObjectsToBeDeleted->Add(line);
618 line->SetLineColor(kBlue);
619 line->SetLineWidth(2);
620 line->SetLineStyle(1);
622 sprintf(c, "Median: %f", median);
623 legend->AddEntry(line, c, "l");
625 for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
626 TLine* linePlusSigma = new TLine(median + nsigma[i] * sigma, 0, median + nsigma[i]*sigma, maxY);
627 //fListOfObjectsToBeDeleted->Add(linePlusSigma);
628 linePlusSigma->SetLineColor(kBlue);
629 linePlusSigma->SetLineStyle(2 + i);
630 linePlusSigma->Draw();
631 TLine* lineMinusSigma = new TLine(median - nsigma[i] * sigma, 0, median - nsigma[i]*sigma, maxY);
632 //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
633 lineMinusSigma->SetLineColor(kBlue);
634 lineMinusSigma->SetLineStyle(2 + i);
635 lineMinusSigma->Draw();
636 sprintf(c, "%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i] * sigma));
637 legend->AddEntry(lineMinusSigma, c, "l");
643 Double_t ltm = GetLTM(entries, values, <mRms, ltmFraction);
644 TLine* line = new TLine(ltm, 0, ltm, maxY);
645 //fListOfObjectsToBeDeleted->Add(line);
646 line->SetLineColor(kGreen+2);
647 line->SetLineWidth(2);
648 line->SetLineStyle(1);
650 sprintf(c, "LTM: %f", ltm);
651 legend->AddEntry(line, c, "l");
653 for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
654 TLine* linePlusSigma = new TLine(ltm + nsigma[i] * ltmRms, 0, ltm + nsigma[i] * ltmRms, maxY);
655 //fListOfObjectsToBeDeleted->Add(linePlusSigma);
656 linePlusSigma->SetLineColor(kGreen+2);
657 linePlusSigma->SetLineStyle(2+i);
658 linePlusSigma->Draw();
660 TLine* lineMinusSigma = new TLine(ltm - nsigma[i] * ltmRms, 0, ltm - nsigma[i] * ltmRms, maxY);
661 //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
662 lineMinusSigma->SetLineColor(kGreen+2);
663 lineMinusSigma->SetLineStyle(2+i);
664 lineMinusSigma->Draw();
665 sprintf(c, "%i #sigma = %f", (Int_t)(nsigma[i]), (Float_t)(nsigma[i] * ltmRms));
666 legend->AddEntry(lineMinusSigma, c, "l");
669 if (!plotMean && !plotMedian && !plotLTM) return -1;
671 gStyle->SetOptStat(oldOptStat);
676 Int_t AliTPCCalibViewer::SigmaCut(const char* drawCommand, Int_t sector, const char* cuts, Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, Bool_t pm, const char *sigmas, Float_t sigmaStep) const {
678 // 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
679 // The data of the distribution Begin_Latex f(x, #mu, #sigma) End_Latex are given in 'array', 'n' specifies the length of the array
680 // '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
681 // 'nbins': number of bins, 'binLow': first bin, 'binUp': last bin
682 // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, Begin_Latex t #sigma End_Latex)
683 // sigmaStep: the binsize of the generated histogram
685 // f(x, #mu, #sigma) #Rightarrow S(t, #mu, #sigma) = #frac{#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 }
689 // Creates a histogram, where you can see, how much of the data are inside sigma-intervals
690 // around the mean/median/LTM
691 // with drawCommand, sector and cuts you specify your input data, see EasyDraw
692 // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma)
693 // sigmaStep: the binsize of the generated histogram
694 // plotMean/plotMedian/plotLTM: specifies where to put the center
696 if (sector >= 0 && sector < 72) {
698 sprintf(sectorChr, "%i", sector);
699 return SigmaCut(drawCommand, sectorChr, cuts, sigmaMax, plotMean, plotMedian, plotLTM, pm, sigmas, sigmaStep);
701 Error("SigmaCut","The TPC contains only sectors between 0 and 71.");
706 Int_t AliTPCCalibViewer::SigmaCut(const char* drawCommand, const char* sector, const char* cuts, Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, Bool_t pm, const char *sigmas, Float_t sigmaStep) const {
708 // Creates a histogram, where you can see, how much of the data are inside sigma-intervals
709 // around the mean/median/LTM
710 // with drawCommand, sector and cuts you specify your input data, see EasyDraw
711 // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma)
712 // sigmaStep: the binsize of the generated histogram
713 // plotMean/plotMedian/plotLTM: specifies where to put the center
716 Double_t ltmFraction = 0.8;
718 TString drawStr(drawCommand);
719 Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
720 if (dangerousToDraw) {
721 Warning("SigmaCut", "The draw string must not contain ':' or '>>'.");
724 drawStr += " >> tempHist";
726 Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
727 TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
728 // FIXME is this histogram deleted automatically?
729 Double_t *values = fTree->GetV1(); // value is the array containing 'entries' numbers
731 Double_t mean = TMath::Mean(entries, values);
732 Double_t median = TMath::Median(entries, values);
733 Double_t sigma = TMath::RMS(entries, values);
735 TLegend * legend = new TLegend(.7,.7, .99, .99, "Cumulative");
736 //fListOfObjectsToBeDeleted->Add(legend);
737 TH1F *cutHistoMean = 0;
738 TH1F *cutHistoMedian = 0;
739 TH1F *cutHistoLTM = 0;
741 TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
742 TVectorF nsigma(sigmasTokens->GetEntriesFast());
743 for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
744 TString str(((TObjString*)sigmasTokens->At(i))->GetString());
745 Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
750 cutHistoMean = AliTPCCalibViewer::SigmaCut(htemp, mean, sigma, sigmaMax, sigmaStep, pm);
752 //fListOfObjectsToBeDeleted->Add(cutHistoMean);
753 cutHistoMean->SetLineColor(kRed);
754 legend->AddEntry(cutHistoMean, "Mean", "l");
755 cutHistoMean->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
756 cutHistoMean->Draw();
757 DrawLines(cutHistoMean, nsigma, legend, kRed, pm);
758 } // if (cutHistoMean)
762 cutHistoMedian = AliTPCCalibViewer::SigmaCut(htemp, median, sigma, sigmaMax, sigmaStep, pm);
763 if (cutHistoMedian) {
764 //fListOfObjectsToBeDeleted->Add(cutHistoMedian);
765 cutHistoMedian->SetLineColor(kBlue);
766 legend->AddEntry(cutHistoMedian, "Median", "l");
767 cutHistoMedian->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
768 if (plotMean && cutHistoMean) cutHistoMedian->Draw("same");
769 else cutHistoMedian->Draw();
770 DrawLines(cutHistoMedian, nsigma, legend, kBlue, pm);
771 } // if (cutHistoMedian)
775 Double_t ltm = GetLTM(entries, values, <mRms, ltmFraction);
776 cutHistoLTM = AliTPCCalibViewer::SigmaCut(htemp, ltm, ltmRms, sigmaMax, sigmaStep, pm);
778 //fListOfObjectsToBeDeleted->Add(cutHistoLTM);
779 cutHistoLTM->SetLineColor(kGreen+2);
780 legend->AddEntry(cutHistoLTM, "LTM", "l");
781 cutHistoLTM->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
782 if (plotMean && cutHistoMean || plotMedian && cutHistoMedian) cutHistoLTM->Draw("same");
783 else cutHistoLTM->Draw();
784 DrawLines(cutHistoLTM, nsigma, legend, kGreen+2, pm);
787 if (!plotMean && !plotMedian && !plotLTM) return -1;
793 Int_t AliTPCCalibViewer::SigmaCutNew(const char* drawCommand, const char* sector, const char* cuts, Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, Bool_t pm, const char *sigmas, Float_t sigmaStep) const {
795 // Creates a histogram, where you can see, how much of the data are inside sigma-intervals
796 // around the mean/median/LTM
797 // with drawCommand, sector and cuts you specify your input data, see EasyDraw
798 // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma)
799 // sigmaStep: the binsize of the generated histogram
800 // plotMean/plotMedian/plotLTM: specifies where to put the center
803 // Double_t ltmFraction = 0.8; //unused
804 // avoid compiler warnings:
807 sigmaStep = sigmaStep;
809 TString drawStr(drawCommand);
810 drawStr += " >> tempHist";
812 Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
813 TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
814 TGraph *cutGraphMean = 0;
815 // TGraph *cutGraphMedian = 0;
816 // TGraph *cutGraphLTM = 0;
817 Double_t *values = fTree->GetV1(); // value is the array containing 'entries' numbers
818 Int_t *index = new Int_t[entries];
819 Float_t *xarray = new Float_t[entries];
820 Float_t *yarray = new Float_t[entries];
821 TMath::Sort(entries, values, index, kFALSE);
823 Double_t mean = TMath::Mean(entries, values);
824 // Double_t median = TMath::Median(entries, values);
825 Double_t sigma = TMath::RMS(entries, values);
827 TLegend * legend = new TLegend(.7,.7, .99, .99, "Cumulative");
828 //fListOfObjectsToBeDeleted->Add(legend);
830 // parse sigmas string
831 TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
832 TVectorF nsigma(sigmasTokens->GetEntriesFast());
833 for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
834 TString str(((TObjString*)sigmasTokens->At(i))->GetString());
835 Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
840 for (Int_t i = 0; i < entries; i++) {
841 xarray[i] = TMath::Abs(values[index[i]] - mean) / sigma;
842 yarray[i] = float(i) / float(entries);
844 cutGraphMean = new TGraph(entries, xarray, yarray);
846 //fListOfObjectsToBeDeleted->Add(cutGraphMean);
847 cutGraphMean->SetLineColor(kRed);
848 legend->AddEntry(cutGraphMean, "Mean", "l");
849 cutGraphMean->SetTitle(Form("%s, Cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
850 cutGraphMean->Draw("alu");
851 DrawLines(cutGraphMean, nsigma, legend, kRed, kTRUE);
856 cutHistoMedian = AliTPCCalibViewer::SigmaCut(htemp, median, sigma, sigmaMax, sigmaStep, pm);
857 if (cutHistoMedian) {
858 fListOfObjectsToBeDeleted->Add(cutHistoMedian);
859 cutHistoMedian->SetLineColor(kBlue);
860 legend->AddEntry(cutHistoMedian, "Median", "l");
861 cutHistoMedian->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
862 if (plotMean && cutHistoMean) cutHistoMedian->Draw("same");
863 else cutHistoMedian->Draw();
864 DrawLines(cutHistoMedian, nsigma, legend, kBlue, pm);
865 } // if (cutHistoMedian)
869 Double_t ltm = GetLTM(entries, values, <mRms, ltmFraction);
870 cutHistoLTM = AliTPCCalibViewer::SigmaCut(htemp, ltm, ltmRms, sigmaMax, sigmaStep, pm);
872 fListOfObjectsToBeDeleted->Add(cutHistoLTM);
873 cutHistoLTM->SetLineColor(kGreen+2);
874 legend->AddEntry(cutHistoLTM, "LTM", "l");
875 cutHistoLTM->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
876 if (plotMean && cutHistoMean || plotMedian && cutHistoMedian) cutHistoLTM->Draw("same");
877 else cutHistoLTM->Draw();
878 DrawLines(cutHistoLTM, nsigma, legend, kGreen+2, pm);
881 if (!plotMean && !plotMedian && !plotLTM) return -1;
887 Int_t AliTPCCalibViewer::Integrate(const char* drawCommand, Int_t sector, const char* cuts, Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, const char *sigmas, Float_t sigmaStep) const {
889 // 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"
890 // "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
891 // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate
892 // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
893 // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
894 // The actual work is done on the array.
896 f(x, #mu, #sigma) #Rightarrow S(t, #mu, #sigma) = #frac{#int_{-#infty}^{#mu + t #sigma} f(x, #mu, #sigma) dx}{ #int_{-#infty}^{+#infty} f(x, #mu, #sigma) dx }
899 if (sector >= 0 && sector < 72) {
901 sprintf(sectorChr, "%i", sector);
902 return Integrate(drawCommand, sectorChr, cuts, sigmaMax, plotMean, plotMedian, plotLTM, sigmas, sigmaStep);
904 Error("Integrate","The TPC contains only sectors between 0 and 71.");
910 Int_t AliTPCCalibViewer::IntegrateOld(const char* drawCommand, const char* sector, const char* cuts, Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, const char *sigmas, Float_t sigmaStep) const {
912 // 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"
913 // "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
914 // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate
915 // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
916 // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
917 // The actual work is done on the array.
919 f(x, #mu, #sigma) #Rightarrow S(t, #mu, #sigma) = #frac{#int_{-#infty}^{#mu + t #sigma} f(x, #mu, #sigma) dx}{ #int_{-#infty}^{+#infty} f(x, #mu, #sigma) dx }
923 Double_t ltmFraction = 0.8;
925 TString drawStr(drawCommand);
926 drawStr += " >> tempHist";
928 Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
929 TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
930 // FIXME is this histogram deleted automatically?
931 Double_t *values = fTree->GetV1(); // value is the array containing 'entries' numbers
933 Double_t mean = TMath::Mean(entries, values);
934 Double_t median = TMath::Median(entries, values);
935 Double_t sigma = TMath::RMS(entries, values);
937 TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
938 TVectorF nsigma(sigmasTokens->GetEntriesFast());
939 for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
940 TString str(((TObjString*)sigmasTokens->At(i))->GetString());
941 Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
945 TLegend * legend = new TLegend(.7,.7, .99, .99, "Integrated histogram");
946 //fListOfObjectsToBeDeleted->Add(legend);
947 TH1F *integralHistoMean = 0;
948 TH1F *integralHistoMedian = 0;
949 TH1F *integralHistoLTM = 0;
952 integralHistoMean = AliTPCCalibViewer::Integrate(htemp, mean, sigma, sigmaMax, sigmaStep);
953 if (integralHistoMean) {
954 //fListOfObjectsToBeDeleted->Add(integralHistoMean);
955 integralHistoMean->SetLineColor(kRed);
956 legend->AddEntry(integralHistoMean, "Mean", "l");
957 integralHistoMean->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
958 integralHistoMean->Draw();
959 DrawLines(integralHistoMean, nsigma, legend, kRed, kTRUE);
963 integralHistoMedian = AliTPCCalibViewer::Integrate(htemp, median, sigma, sigmaMax, sigmaStep);
964 if (integralHistoMedian) {
965 //fListOfObjectsToBeDeleted->Add(integralHistoMedian);
966 integralHistoMedian->SetLineColor(kBlue);
967 legend->AddEntry(integralHistoMedian, "Median", "l");
968 integralHistoMedian->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
969 if (plotMean && integralHistoMean) integralHistoMedian->Draw("same");
970 else integralHistoMedian->Draw();
971 DrawLines(integralHistoMedian, nsigma, legend, kBlue, kTRUE);
976 Double_t ltm = GetLTM(entries, values, <mRms, ltmFraction);
977 integralHistoLTM = AliTPCCalibViewer::Integrate(htemp, ltm, ltmRms, sigmaMax, sigmaStep);
978 if (integralHistoLTM) {
979 //fListOfObjectsToBeDeleted->Add(integralHistoLTM);
980 integralHistoLTM->SetLineColor(kGreen+2);
981 legend->AddEntry(integralHistoLTM, "LTM", "l");
982 integralHistoLTM->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
983 if (plotMean && integralHistoMean || plotMedian && integralHistoMedian) integralHistoLTM->Draw("same");
984 else integralHistoLTM->Draw();
985 DrawLines(integralHistoLTM, nsigma, legend, kGreen+2, kTRUE);
988 if (!plotMean && !plotMedian && !plotLTM) return -1;
994 Int_t AliTPCCalibViewer::Integrate(const char* drawCommand, const char* sector, const char* cuts, Float_t sigmaMax, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM, const char *sigmas, Float_t sigmaStep) const {
996 // 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"
997 // "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
998 // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate
999 // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
1000 // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
1001 // The actual work is done on the array.
1003 f(x, #mu, #sigma) #Rightarrow S(t, #mu, #sigma) = #frac{#int_{-#infty}^{#mu + t #sigma} f(x, #mu, #sigma) dx}{ #int_{-#infty}^{+#infty} f(x, #mu, #sigma) dx }
1007 Double_t ltmFraction = 0.8;
1008 // avoid compiler warnings:
1009 sigmaMax = sigmaMax;
1010 sigmaStep = sigmaStep;
1012 TString drawStr(drawCommand);
1013 Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
1014 if (dangerousToDraw) {
1015 Warning("Integrate", "The draw string must not contain ':' or '>>'.");
1018 drawStr += " >> tempHist";
1020 Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
1021 TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
1022 TGraph *integralGraphMean = 0;
1023 TGraph *integralGraphMedian = 0;
1024 TGraph *integralGraphLTM = 0;
1025 Double_t *values = fTree->GetV1(); // value is the array containing 'entries' numbers
1026 Int_t *index = new Int_t[entries];
1027 Float_t *xarray = new Float_t[entries];
1028 Float_t *yarray = new Float_t[entries];
1029 TMath::Sort(entries, values, index, kFALSE);
1031 Double_t mean = TMath::Mean(entries, values);
1032 Double_t median = TMath::Median(entries, values);
1033 Double_t sigma = TMath::RMS(entries, values);
1035 // parse sigmas string
1036 TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");
1037 TVectorF nsigma(sigmasTokens->GetEntriesFast());
1038 for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
1039 TString str(((TObjString*)sigmasTokens->At(i))->GetString());
1040 Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
1044 TLegend * legend = new TLegend(.7,.7, .99, .99, "Integrated histogram");
1045 //fListOfObjectsToBeDeleted->Add(legend);
1048 for (Int_t i = 0; i < entries; i++) {
1049 xarray[i] = (values[index[i]] - mean) / sigma;
1050 yarray[i] = float(i) / float(entries);
1052 integralGraphMean = new TGraph(entries, xarray, yarray);
1053 if (integralGraphMean) {
1054 //fListOfObjectsToBeDeleted->Add(integralGraphMean);
1055 integralGraphMean->SetLineColor(kRed);
1056 legend->AddEntry(integralGraphMean, "Mean", "l");
1057 integralGraphMean->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
1058 integralGraphMean->Draw("alu");
1059 DrawLines(integralGraphMean, nsigma, legend, kRed, kTRUE);
1063 for (Int_t i = 0; i < entries; i++) {
1064 xarray[i] = (values[index[i]] - median) / sigma;
1065 yarray[i] = float(i) / float(entries);
1067 integralGraphMedian = new TGraph(entries, xarray, yarray);
1068 if (integralGraphMedian) {
1069 //fListOfObjectsToBeDeleted->Add(integralGraphMedian);
1070 integralGraphMedian->SetLineColor(kBlue);
1071 legend->AddEntry(integralGraphMedian, "Median", "l");
1072 integralGraphMedian->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
1073 if (plotMean && integralGraphMean) integralGraphMedian->Draw("samelu");
1074 else integralGraphMedian->Draw("alu");
1075 DrawLines(integralGraphMedian, nsigma, legend, kBlue, kTRUE);
1079 Double_t ltmRms = 0;
1080 Double_t ltm = GetLTM(entries, values, <mRms, ltmFraction);
1081 for (Int_t i = 0; i < entries; i++) {
1082 xarray[i] = (values[index[i]] - ltm) / ltmRms;
1083 yarray[i] = float(i) / float(entries);
1085 integralGraphLTM = new TGraph(entries, xarray, yarray);
1086 if (integralGraphLTM) {
1087 //fListOfObjectsToBeDeleted->Add(integralGraphLTM);
1088 integralGraphLTM->SetLineColor(kGreen+2);
1089 legend->AddEntry(integralGraphLTM, "LTM", "l");
1090 integralGraphLTM->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
1091 if (plotMean && integralGraphMean || plotMedian && integralGraphMedian) integralGraphLTM->Draw("samelu");
1092 else integralGraphLTM->Draw("alu");
1093 DrawLines(integralGraphLTM, nsigma, legend, kGreen+2, kTRUE);
1096 if (!plotMean && !plotMedian && !plotLTM) return -1;
1102 void AliTPCCalibViewer::DrawLines(TH1F *histogram, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const {
1104 // Private function for SigmaCut(...) and Integrate(...)
1105 // Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
1108 // start to draw the lines, loop over requested sigmas
1110 for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
1112 Int_t bin = histogram->GetXaxis()->FindBin(nsigma[i]);
1113 TLine* lineUp = new TLine(nsigma[i], 0, nsigma[i], histogram->GetBinContent(bin));
1114 //fListOfObjectsToBeDeleted->Add(lineUp);
1115 lineUp->SetLineColor(color);
1116 lineUp->SetLineStyle(2 + i);
1118 TLine* lineLeft = new TLine(nsigma[i], histogram->GetBinContent(bin), 0, histogram->GetBinContent(bin));
1119 //fListOfObjectsToBeDeleted->Add(lineLeft);
1120 lineLeft->SetLineColor(color);
1121 lineLeft->SetLineStyle(2 + i);
1123 sprintf(c, "Fraction(%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin));
1124 legend->AddEntry(lineLeft, c, "l");
1127 Int_t bin = histogram->GetXaxis()->FindBin(nsigma[i]);
1128 TLine* lineUp1 = new TLine(nsigma[i], 0, nsigma[i], histogram->GetBinContent(bin));
1129 //fListOfObjectsToBeDeleted->Add(lineUp1);
1130 lineUp1->SetLineColor(color);
1131 lineUp1->SetLineStyle(2 + i);
1133 TLine* lineLeft1 = new TLine(nsigma[i], histogram->GetBinContent(bin), histogram->GetBinLowEdge(0)+histogram->GetBinWidth(0), histogram->GetBinContent(bin));
1134 //fListOfObjectsToBeDeleted->Add(lineLeft1);
1135 lineLeft1->SetLineColor(color);
1136 lineLeft1->SetLineStyle(2 + i);
1138 sprintf(c, "Fraction(+%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin));
1139 legend->AddEntry(lineLeft1, c, "l");
1140 bin = histogram->GetXaxis()->FindBin(-nsigma[i]);
1141 TLine* lineUp2 = new TLine(-nsigma[i], 0, -nsigma[i], histogram->GetBinContent(bin));
1142 //fListOfObjectsToBeDeleted->Add(lineUp2);
1143 lineUp2->SetLineColor(color);
1144 lineUp2->SetLineStyle(2 + i);
1146 TLine* lineLeft2 = new TLine(-nsigma[i], histogram->GetBinContent(bin), histogram->GetBinLowEdge(0)+histogram->GetBinWidth(0), histogram->GetBinContent(bin));
1147 //fListOfObjectsToBeDeleted->Add(lineLeft2);
1148 lineLeft2->SetLineColor(color);
1149 lineLeft2->SetLineStyle(2 + i);
1151 sprintf(c, "Fraction(-%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin));
1152 legend->AddEntry(lineLeft2, c, "l");
1154 } // for (Int_t i = 0; i < nsigma.GetNoElements(); i++)
1158 void AliTPCCalibViewer::DrawLines(TGraph *graph, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const {
1160 // Private function for SigmaCut(...) and Integrate(...)
1161 // Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
1164 // start to draw the lines, loop over requested sigmas
1166 for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
1168 TLine* lineUp = new TLine(nsigma[i], 0, nsigma[i], graph->Eval(nsigma[i]));
1169 //fListOfObjectsToBeDeleted->Add(lineUp);
1170 lineUp->SetLineColor(color);
1171 lineUp->SetLineStyle(2 + i);
1173 TLine* lineLeft = new TLine(nsigma[i], graph->Eval(nsigma[i]), 0, graph->Eval(nsigma[i]));
1174 //fListOfObjectsToBeDeleted->Add(lineLeft);
1175 lineLeft->SetLineColor(color);
1176 lineLeft->SetLineStyle(2 + i);
1178 sprintf(c, "Fraction(%f #sigma) = %f",nsigma[i], graph->Eval(nsigma[i]));
1179 legend->AddEntry(lineLeft, c, "l");
1182 TLine* lineUp1 = new TLine(nsigma[i], 0, nsigma[i], graph->Eval(nsigma[i]));
1183 //fListOfObjectsToBeDeleted->Add(lineUp1);
1184 lineUp1->SetLineColor(color);
1185 lineUp1->SetLineStyle(2 + i);
1187 TLine* lineLeft1 = new TLine(nsigma[i], graph->Eval(nsigma[i]), graph->GetHistogram()->GetXaxis()->GetBinLowEdge(0), graph->Eval(nsigma[i]));
1188 //fListOfObjectsToBeDeleted->Add(lineLeft1);
1189 lineLeft1->SetLineColor(color);
1190 lineLeft1->SetLineStyle(2 + i);
1192 sprintf(c, "Fraction(+%f #sigma) = %f",nsigma[i], graph->Eval(nsigma[i]));
1193 legend->AddEntry(lineLeft1, c, "l");
1194 TLine* lineUp2 = new TLine(-nsigma[i], 0, -nsigma[i], graph->Eval(-nsigma[i]));
1195 //fListOfObjectsToBeDeleted->Add(lineUp2);
1196 lineUp2->SetLineColor(color);
1197 lineUp2->SetLineStyle(2 + i);
1199 TLine* lineLeft2 = new TLine(-nsigma[i], graph->Eval(-nsigma[i]), graph->GetHistogram()->GetXaxis()->GetBinLowEdge(0), graph->Eval(-nsigma[i]));
1200 //fListOfObjectsToBeDeleted->Add(lineLeft2);
1201 lineLeft2->SetLineColor(color);
1202 lineLeft2->SetLineStyle(2 + i);
1204 sprintf(c, "Fraction(-%f #sigma) = %f",nsigma[i], graph->Eval(-nsigma[i]));
1205 legend->AddEntry(lineLeft2, c, "l");
1207 } // for (Int_t i = 0; i < nsigma.GetNoElements(); i++)
1219 Int_t AliTPCCalibViewer::GetBin(Float_t value, Int_t nbins, Double_t binLow, Double_t binUp){
1220 // Returns the 'bin' for 'value'
1221 // The interval between 'binLow' and 'binUp' is divided into 'nbins' equidistant bins
1222 // avoid index out of bounds error: 'if (bin < binLow) bin = binLow' and vice versa
1224 GetBin(value) = #frac{nbins - 1}{binUp - binLow} #upoint (value - binLow) +1
1228 Int_t bin = TMath::Nint( (Float_t)(value - binLow) / (Float_t)(binUp - binLow) * (nbins-1) ) + 1;
1229 // avoid index out of bounds:
1230 if (value < binLow) bin = 0;
1231 if (value > binUp) bin = nbins + 1;
1237 Double_t AliTPCCalibViewer::GetLTM(Int_t n, Double_t *array, Double_t *sigma, Double_t fraction){
1239 // returns the LTM and sigma
1241 Double_t *ddata = new Double_t[n];
1242 Double_t mean = 0, lsigma = 0;
1244 for (UInt_t i = 0; i < (UInt_t)n; i++) {
1245 ddata[nPoints]= array[nPoints];
1248 Int_t hh = TMath::Min(TMath::Nint(fraction * nPoints), Int_t(n));
1249 AliMathBase::EvaluateUni(nPoints, ddata, mean, lsigma, hh);
1250 if (sigma) *sigma = lsigma;
1256 TH1F* AliTPCCalibViewer::SigmaCut(TH1F *histogram, Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep, Bool_t pm) {
1258 // 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
1259 // The data of the distribution Begin_Latex f(x, #mu, #sigma) End_Latex are given in 'histogram'
1260 // '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
1261 // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, Begin_Latex t #sigma End_Latex)
1262 // sigmaStep: the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
1263 // pm: Decide weather Begin_Latex t > 0 End_Latex (first case) or Begin_Latex t End_Latex arbitrary (secound case)
1264 // The actual work is done on the array.
1266 f(x, #mu, #sigma) #Rightarrow S(t, #mu, #sigma) = #frac{#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
1268 f(x, #mu, #sigma) #Rightarrow S(t, #mu, #sigma) = #frac{#int_{#mu}^{#mu + t #sigma} f(x, #mu, #sigma) dx}{ #int_{-#infty}^{+#infty} f(x, #mu, #sigma) dx }
1273 Float_t sigma = 1.5;
1274 Float_t sigmaMax = 4;
1275 gROOT->SetStyle("Plain");
1276 TH1F *distribution = new TH1F("Distribution1", "Distribution f(x, #mu, #sigma)", 1000,-5,5);
1278 for (Int_t i = 0; i <50000;i++) distribution->Fill(rand.Gaus(mean, sigma));
1279 Float_t *ar = distribution->GetArray();
1281 TCanvas* macro_example_canvas = new TCanvas("macro_example_canvas_SigmaCut", "", 350, 350);
1282 macro_example_canvas->Divide(0,3);
1283 TVirtualPad *pad1 = macro_example_canvas->cd(1);
1286 distribution->Draw();
1287 TVirtualPad *pad2 = macro_example_canvas->cd(2);
1291 TH1F *shist = AliTPCCalibViewer::SigmaCut(distribution, mean, sigma, sigmaMax);
1292 shist->SetNameTitle("Cumulative","Cumulative S(t, #mu, #sigma)");
1294 TVirtualPad *pad3 = macro_example_canvas->cd(3);
1297 TH1F *shistPM = AliTPCCalibViewer::SigmaCut(distribution, mean, sigma, sigmaMax, -1, kTRUE);
1299 return macro_example_canvas;
1304 Float_t *array = histogram->GetArray();
1305 Int_t nbins = histogram->GetXaxis()->GetNbins();
1306 Float_t binLow = histogram->GetXaxis()->GetXmin();
1307 Float_t binUp = histogram->GetXaxis()->GetXmax();
1308 return AliTPCCalibViewer::SigmaCut(nbins, array, mean, sigma, nbins, binLow, binUp, sigmaMax, sigmaStep, pm);
1312 TH1F* AliTPCCalibViewer::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){
1314 // 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
1315 // The data of the distribution Begin_Latex f(x, #mu, #sigma) End_Latex are given in 'array', 'n' specifies the length of the array
1316 // '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
1317 // 'nbins': number of bins, 'binLow': first bin, 'binUp': last bin
1318 // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, Begin_Latex t #sigma End_Latex)
1319 // sigmaStep: the binsize of the generated histogram
1320 // Here the actual work is done.
1322 if (sigma == 0) return 0;
1323 Float_t binWidth = (binUp-binLow)/(nbins - 1);
1324 if (sigmaStep <= 0) sigmaStep = binWidth;
1325 Int_t kbins = (Int_t)(sigmaMax * sigma / sigmaStep) + 1; // + 1 due to overflow bin in histograms
1326 if (pm) kbins = 2 * (Int_t)(sigmaMax * sigma / sigmaStep) + 1;
1327 Float_t kbinLow = !pm ? 0 : -sigmaMax;
1328 Float_t kbinUp = sigmaMax;
1329 TH1F *hist = new TH1F("sigmaCutHisto","Cumulative; Multiples of #sigma; Fraction of included data", kbins, kbinLow, kbinUp);
1330 hist->SetDirectory(0);
1333 // calculate normalization
1334 Double_t normalization = 0;
1335 for (Int_t i = 0; i <= n; i++) {
1336 normalization += array[i];
1339 // given units: units from given histogram
1340 // sigma units: in units of sigma
1341 // iDelta: integrate in interval (mean +- iDelta), given units
1342 // x: ofset from mean for integration, given units
1345 // printf("nbins: %i, binLow: %f, binUp: %f \n", nbins, binLow, binUp);
1347 for (Float_t iDelta = 0; iDelta <= sigmaMax * sigma; iDelta += sigmaStep) {
1349 Double_t valueP = array[GetBin(mean, nbins, binLow, binUp)];
1350 Double_t valueM = array[GetBin(mean-binWidth, nbins, binLow, binUp)];
1351 // add bin of mean value only once to the histogram
1352 // printf("++ adding bins: ");
1353 for (Float_t x = binWidth; x <= iDelta; x += binWidth) {
1354 valueP += (mean + x <= binUp) ? array[GetBin(mean + x, nbins, binLow, binUp)] : 0;
1355 valueM += (mean-binWidth - x >= binLow) ? array[GetBin(mean-binWidth - x, nbins, binLow, binUp)] : 0;
1356 // printf("%i, ", GetBin(mean + x, nbins, binLow, binUp));
1359 if (valueP / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail +++ \n", valueP, normalization);
1360 if (valueP / normalization > 100) return hist;
1361 if (valueM / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail +++ \n", valueM, normalization);
1362 if (valueM / normalization > 100) return hist;
1363 valueP = (valueP / normalization);
1364 valueM = (valueM / normalization);
1366 Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
1367 hist->SetBinContent(bin, valueP);
1368 bin = GetBin(-iDelta/sigma, kbins, kbinLow, kbinUp);
1369 hist->SetBinContent(bin, valueM);
1372 Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
1373 hist->SetBinContent(bin, valueP + valueM);
1374 // printf(" first integration bin: %i, last integration bin in + direction: %i \n", GetBin(mean+binWidth, nbins, binLow, binUp), GetBin(iDelta, nbins, binLow, binUp));
1375 // printf(" first integration bin: %i, last integration bin in - direction: %i \n", GetBin(mean+binWidth, nbins, binLow, binUp), GetBin(-iDelta, nbins, binLow, binUp));
1376 // printf(" value: %f, normalization: %f, iDelta: %f, Bin: %i \n", valueP+valueM, normalization, iDelta, bin);
1379 //hist->SetMaximum(0.7);
1380 if (!pm) hist->SetMaximum(1.2);
1385 TH1F* AliTPCCalibViewer::SigmaCut(Int_t n, Double_t *array, Double_t mean, Double_t sigma, Int_t nbins, Double_t *xbins, Double_t sigmaMax){
1387 // SigmaCut for variable binsize
1388 // NOT YET IMPLEMENTED !!!
1390 printf("SigmaCut with variable binsize, Not yet implemented\n");
1391 // avoid compiler warnings:
1404 TH1F* AliTPCCalibViewer::Integrate(TH1F *histogram, Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep){
1406 // 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"
1407 // "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
1408 // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate
1409 // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
1410 // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
1411 // The actual work is done on the array.
1413 f(x, #mu, #sigma) #Rightarrow S(t, #mu, #sigma) = #frac{#int_{-#infty}^{#mu + t #sigma} f(x, #mu, #sigma) dx}{ #int_{-#infty}^{+#infty} f(x, #mu, #sigma) dx }
1418 Float_t sigma = 1.5;
1419 Float_t sigmaMax = 4;
1420 gROOT->SetStyle("Plain");
1421 TH1F *distribution = new TH1F("Distribution2", "Distribution f(x, #mu, #sigma)", 1000,-5,5);
1423 for (Int_t i = 0; i <50000;i++) distribution->Fill(rand.Gaus(mean, sigma));
1424 Float_t *ar = distribution->GetArray();
1426 TCanvas* macro_example_canvas = new TCanvas("macro_example_canvas_Integrate", "", 350, 350);
1427 macro_example_canvas->Divide(0,2);
1428 TVirtualPad *pad1 = macro_example_canvas->cd(1);
1431 distribution->Draw();
1432 TVirtualPad *pad2 = macro_example_canvas->cd(2);
1435 TH1F *shist = AliTPCCalibViewer::Integrate(distribution, mean, sigma, sigmaMax);
1436 shist->SetNameTitle("Cumulative","Cumulative S(t, #mu, #sigma)");
1439 return macro_example_canvas_Integrate;
1445 Float_t *array = histogram->GetArray();
1446 Int_t nbins = histogram->GetXaxis()->GetNbins();
1447 Float_t binLow = histogram->GetXaxis()->GetXmin();
1448 Float_t binUp = histogram->GetXaxis()->GetXmax();
1449 return AliTPCCalibViewer::Integrate(nbins, array, nbins, binLow, binUp, mean, sigma, sigmaMax, sigmaStep);
1453 TH1F* AliTPCCalibViewer::Integrate(Int_t n, Float_t *array, Int_t nbins, Float_t binLow, Float_t binUp, Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep){
1454 // 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"
1455 // "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
1456 // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate
1457 // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
1458 // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
1459 // Here the actual work is done.
1461 Bool_t givenUnits = kTRUE;
1462 if (sigma != 0 && sigmaMax != 0) givenUnits = kFALSE;
1465 sigmaMax = (binUp - binLow) / 2.;
1468 Float_t binWidth = (binUp-binLow)/(nbins - 1);
1469 if (sigmaStep <= 0) sigmaStep = binWidth;
1470 Int_t kbins = (Int_t)(sigmaMax * sigma / sigmaStep) + 1; // + 1 due to overflow bin in histograms
1471 Float_t kbinLow = givenUnits ? binLow : -sigmaMax;
1472 Float_t kbinUp = givenUnits ? binUp : sigmaMax;
1474 if (givenUnits) hist = new TH1F("integratedHisto","Integrated Histogram; Given x; Fraction of included data", kbins, kbinLow, kbinUp);
1475 if (!givenUnits) hist = new TH1F("integratedHisto","Integrated Histogram; Multiples of #sigma; Fraction of included data", kbins, kbinLow, kbinUp);
1476 hist->SetDirectory(0);
1479 // calculate normalization
1480 // printf("calculating normalization, integrating from bin 1 to %i \n", n);
1481 Double_t normalization = 0;
1482 for (Int_t i = 1; i <= n; i++) {
1483 normalization += array[i];
1485 // printf("normalization: %f \n", normalization);
1487 // given units: units from given histogram
1488 // sigma units: in units of sigma
1489 // iDelta: integrate in interval (mean +- iDelta), given units
1490 // x: ofset from mean for integration, given units
1494 for (Float_t iDelta = mean - sigmaMax * sigma; iDelta <= mean + sigmaMax * sigma; iDelta += sigmaStep) {
1497 for (Float_t x = mean - sigmaMax * sigma; x <= iDelta; x += binWidth) {
1498 value += (x <= binUp && x >= binLow) ? array[GetBin(x, nbins, binLow, binUp)] : 0;
1500 if (value / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail +++ \n", value, normalization);
1501 if (value / normalization > 100) return hist;
1502 Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
1503 // printf("first integration bin: %i, last integration bin: %i \n", GetBin(mean - sigmaMax * sigma, nbins, binLow, binUp), GetBin(iDelta, nbins, binLow, binUp));
1504 // printf("value: %f, normalization: %f, normalized value: %f, iDelta: %f, Bin: %i \n", value, normalization, value/normalization, iDelta, bin);
1505 value = (value / normalization);
1506 hist->SetBinContent(bin, value);
1515 ////////////////////////
1516 // end of Array tools //
1517 ////////////////////////
1521 //_____________________________________________________________________________
1522 AliTPCCalPad* AliTPCCalibViewer::GetCalPad(const char* desiredData, char* cuts, char* calPadName) const {
1524 // creates a AliTPCCalPad out of the 'desiredData'
1525 // the functionality of EasyDraw1D is used
1526 // calPadName specifies the name of the created AliTPCCalPad
1527 // - this takes a while -
1529 TString drawStr(desiredData);
1530 drawStr.Append(":channel");
1531 drawStr.Append(fAbbreviation);
1532 AliTPCCalPad * createdCalPad = new AliTPCCalPad(calPadName, calPadName);
1534 for (Int_t sec = 0; sec < 72; sec++) {
1535 entries = EasyDraw1D(drawStr.Data(), (Int_t)sec, cuts, "goff");
1536 if (entries == -1) return 0;
1537 for (Int_t i = 0; i < entries; i++)
1538 createdCalPad->GetCalROC(sec)->SetValue((UInt_t)(fTree->GetV2()[i]), (Float_t)(fTree->GetV1()[i]));
1540 return createdCalPad;
1543 //_____________________________________________________________________________
1544 AliTPCCalROC* AliTPCCalibViewer::GetCalROC(const char* desiredData, UInt_t sector, char* cuts) const {
1546 // creates a AliTPCCalROC out of the desiredData
1547 // the functionality of EasyDraw1D is used
1548 // sector specifies the sector of the created AliTPCCalROC
1550 TString drawStr(desiredData);
1551 drawStr.Append(":channel");
1552 drawStr.Append(fAbbreviation);
1553 Int_t entries = EasyDraw1D(drawStr.Data(), (Int_t)sector, cuts, "goff");
1554 if (entries == -1) return 0;
1555 AliTPCCalROC * createdROC = new AliTPCCalROC(sector);
1556 for (Int_t i = 0; i < entries; i++)
1557 createdROC->SetValue((UInt_t)(fTree->GetV2()[i]), fTree->GetV1()[i]);
1562 TObjArray* AliTPCCalibViewer::GetListOfVariables(Bool_t printList) {
1564 // scan the tree - produces a list of available variables in the tree
1565 // printList: print the list to the screen, after the scan is done
1567 TObjArray* arr = new TObjArray();
1568 TObjString* str = 0;
1569 Int_t nentries = fTree->GetListOfBranches()->GetEntries();
1570 for (Int_t i = 0; i < nentries; i++) {
1571 str = new TObjString(fTree->GetListOfBranches()->At(i)->GetName());
1572 str->String().ReplaceAll("_Median", "");
1573 str->String().ReplaceAll("_Mean", "");
1574 str->String().ReplaceAll("_RMS", "");
1575 str->String().ReplaceAll("_LTM", "");
1576 str->String().ReplaceAll("_OutlierCutted", "");
1577 str->String().ReplaceAll(".", "");
1578 if (!arr->FindObject(str) &&
1579 !(str->String() == "channel" || str->String() == "gx" || str->String() == "gy" ||
1580 str->String() == "lx" || str->String() == "ly" || str->String() == "pad" ||
1581 str->String() == "row" || str->String() == "rpad" || str->String() == "sector" ))
1585 // loop over all friends (if there are some) and add them to the list
1586 if (fTree->GetListOfFriends()) {
1587 for (Int_t ifriend = 0; ifriend < fTree->GetListOfFriends()->GetEntries(); ifriend++){
1588 // printf("iterating through friendlist, currently at %i\n", ifriend);
1589 // printf("working with %s\n", fTree->GetListOfFriends()->At(ifriend)->ClassName());
1590 if (TString(fTree->GetListOfFriends()->At(ifriend)->ClassName()) != "TFriendElement") continue; // no friendElement found
1591 TFriendElement *friendElement = (TFriendElement*)fTree->GetListOfFriends()->At(ifriend);
1592 if (friendElement->GetTree() == 0) continue; // no tree found in friendElement
1593 // printf("friend found \n");
1594 for (Int_t i = 0; i < friendElement->GetTree()->GetListOfBranches()->GetEntries(); i++) {
1595 // printf("iterating through friendelement entries, currently at %i\n", i);
1596 str = new TObjString(friendElement->GetTree()->GetListOfBranches()->At(i)->GetName());
1597 str->String().ReplaceAll("_Median", "");
1598 str->String().ReplaceAll("_Mean", "");
1599 str->String().ReplaceAll("_RMS", "");
1600 str->String().ReplaceAll("_LTM", "");
1601 str->String().ReplaceAll("_OutlierCutted", "");
1602 str->String().ReplaceAll(".", "");
1603 if (!(str->String() == "channel" || str->String() == "gx" || str->String() == "gy" ||
1604 str->String() == "lx" || str->String() == "ly" || str->String() == "pad" ||
1605 str->String() == "row" || str->String() == "rpad" || str->String() == "sector" )){
1606 // insert "<friendName>." at the beginning: (<friendName> is per default "R")
1607 str->String().Insert(0, ".");
1608 str->String().Insert(0, friendElement->GetName());
1609 if (!arr->FindObject(str)) arr->Add(str);
1610 // printf("added string %s \n", str->String().Data());
1614 } // if (fTree->GetListOfFriends())
1617 // ((TFriendElement*)gui->GetViewer()->GetTree()->GetListOfFriends()->At(0))->GetTree()->GetListOfBranches()->At(0)->GetName()
1618 // ((TFriendElement*)gui->GetViewer()->GetTree()->GetListOfFriends()->At(0))->GetTree()->GetListOfBranches()
1622 TIterator* iter = arr->MakeIterator();
1624 TObjString* currentStr = 0;
1625 while ( (currentStr = (TObjString*)(iter->Next())) ) {
1626 std::cout << currentStr->GetString().Data() << std::endl;
1634 TObjArray* AliTPCCalibViewer::GetListOfNormalizationVariables(Bool_t printList) const{
1636 // produces a list of available variables for normalization in the tree
1637 // printList: print the list to the screen, after the scan is done
1639 TObjArray* arr = new TObjArray();
1640 arr->Add(new TObjString("_Mean"));
1641 arr->Add(new TObjString("_Mean_OutlierCutted"));
1642 arr->Add(new TObjString("_Median"));
1643 arr->Add(new TObjString("_Median_OutlierCutted"));
1644 arr->Add(new TObjString("_LTM"));
1645 arr->Add(new TObjString("_LTM_OutlierCutted"));
1646 arr->Add(new TObjString(Form("LFitIntern_4_8%s", fAppendString.Data())));
1647 arr->Add(new TObjString(Form("GFitIntern_Lin%s", fAppendString.Data())));
1648 arr->Add(new TObjString(Form("GFitIntern_Par%s", fAppendString.Data())));
1649 arr->Add(new TObjString("FitLinLocal"));
1650 arr->Add(new TObjString("FitLinGlobal"));
1651 arr->Add(new TObjString("FitParLocal"));
1652 arr->Add(new TObjString("FitParGlobal"));
1655 TIterator* iter = arr->MakeIterator();
1657 TObjString* currentStr = 0;
1658 while ((currentStr = (TObjString*)(iter->Next()))) {
1659 std::cout << currentStr->GetString().Data() << std::endl;
1667 TFriendElement* AliTPCCalibViewer::AddReferenceTree(const char* filename, const char* treename, const char* refname){
1669 // add a reference tree to the current tree
1670 // by default the treename is 'calPads' and the reference treename is 'R'
1672 TFile *file = new TFile(filename);
1673 fListOfObjectsToBeDeleted->Add(file);
1674 TTree * tree = (TTree*)file->Get(treename);
1675 return AddFriend(tree, refname);
1679 TObjArray* AliTPCCalibViewer::GetArrayOfCalPads(){
1681 // Returns a TObjArray with all AliTPCCalPads that are stored in the tree
1682 // - this takes a while -
1684 TObjArray *listOfCalPads = GetListOfVariables();
1685 TObjArray *calPadsArray = new TObjArray();
1686 Int_t numberOfCalPads = listOfCalPads->GetEntries();
1687 for (Int_t i = 0; i < numberOfCalPads; i++) {
1688 std::cout << "Creating calPad " << (i+1) << " of " << numberOfCalPads << "\r" << std::flush;
1689 char* calPadName = (char*)((TObjString*)(listOfCalPads->At(i)))->GetString().Data();
1690 TString drawCommand = ((TObjString*)(listOfCalPads->At(i)))->GetString();
1691 drawCommand.Append(fAbbreviation.Data());
1692 AliTPCCalPad* calPad = GetCalPad(drawCommand.Data(), "", calPadName);
1693 calPadsArray->Add(calPad);
1695 std::cout << std::endl;
1696 listOfCalPads->Delete();
1697 delete listOfCalPads;
1698 return calPadsArray;
1702 TString* AliTPCCalibViewer::Fit(const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, TVectorD &fitParam, TMatrixD &covMatrix){
1704 // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
1705 // returns chi2, fitParam and covMatrix
1706 // returns TString with fitted formula
1709 TString formulaStr(formula);
1710 TString drawStr(drawCommand);
1711 TString cutStr(cuts);
1714 drawStr.ReplaceAll(fAbbreviation, fAppendString);
1715 cutStr.ReplaceAll(fAbbreviation, fAppendString);
1716 formulaStr.ReplaceAll(fAbbreviation, fAppendString);
1718 formulaStr.ReplaceAll("++", fAbbreviation);
1719 TObjArray* formulaTokens = formulaStr.Tokenize(fAbbreviation.Data());
1720 Int_t dim = formulaTokens->GetEntriesFast();
1722 fitParam.ResizeTo(dim);
1723 covMatrix.ResizeTo(dim,dim);
1725 TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
1726 fitter->StoreData(kTRUE);
1727 fitter->ClearPoints();
1729 Int_t entries = Draw(drawStr.Data(), cutStr.Data(), "goff");
1730 if (entries == -1) return new TString("An ERROR has occured during fitting!");
1731 Double_t **values = new Double_t*[dim+1] ;
1733 for (Int_t i = 0; i < dim + 1; i++){
1735 if (i < dim) centries = fTree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff");
1736 else centries = fTree->Draw(drawStr.Data(), cutStr.Data(), "goff");
1738 if (entries != centries) return new TString("An ERROR has occured during fitting!");
1739 values[i] = new Double_t[entries];
1740 memcpy(values[i], fTree->GetV1(), entries*sizeof(Double_t));
1743 // add points to the fitter
1744 for (Int_t i = 0; i < entries; i++){
1746 for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
1747 fitter->AddPoint(x, values[dim][i], 1);
1751 fitter->GetParameters(fitParam);
1752 fitter->GetCovarianceMatrix(covMatrix);
1753 chi2 = fitter->GetChisquare();
1756 TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula;
1758 for (Int_t iparam = 0; iparam < dim; iparam++) {
1759 returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
1760 if (iparam < dim-1) returnFormula.Append("+");
1762 returnFormula.Append(" )");
1763 delete formulaTokens;
1766 return preturnFormula;
1770 void AliTPCCalibViewer::MakeTreeWithObjects(const char * fileName, TObjArray * array, const char * mapFileName) {
1772 // Write tree with all available information
1773 // im mapFileName is speciefied, the Map information are also written to the tree
1774 // AliTPCCalPad-Objects are written directly to the tree, so that they can be accessd later on
1775 // (does not work!!!)
1777 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
1779 TObjArray* mapIROCs = 0;
1780 TObjArray* mapOROCs = 0;
1781 TVectorF *mapIROCArray = 0;
1782 TVectorF *mapOROCArray = 0;
1783 Int_t mapEntries = 0;
1784 TString* mapNames = 0;
1787 TFile mapFile(mapFileName, "read");
1789 TList* listOfROCs = mapFile.GetListOfKeys();
1790 mapEntries = listOfROCs->GetEntries()/2;
1791 mapIROCs = new TObjArray(mapEntries*2);
1792 mapOROCs = new TObjArray(mapEntries*2);
1793 mapIROCArray = new TVectorF[mapEntries];
1794 mapOROCArray = new TVectorF[mapEntries];
1796 mapNames = new TString[mapEntries];
1797 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1798 TString rocName(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
1799 rocName.Remove(rocName.Length()-4, 4);
1800 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((rocName + "IROC").Data()), ivalue);
1801 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((rocName + "OROC").Data()), ivalue);
1802 mapNames[ivalue].Append(rocName);
1805 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1806 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
1807 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
1809 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
1810 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
1811 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
1812 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
1815 } // if (mapFileName)
1817 TTreeSRedirector cstream(fileName);
1818 Int_t arrayEntries = array->GetEntries();
1820 // Read names of AliTPCCalPads and save them in names[]
1821 TString* names = new TString[arrayEntries];
1822 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
1823 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
1825 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
1827 TVectorF *vectorArray = new TVectorF[arrayEntries];
1828 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
1829 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
1833 // fill vectors of variable per pad
1835 TVectorF *posArray = new TVectorF[8];
1836 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
1837 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
1839 Float_t posG[3] = {0};
1840 Float_t posL[3] = {0};
1842 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
1843 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
1844 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
1845 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
1846 posArray[0][ichannel] = irow;
1847 posArray[1][ichannel] = ipad;
1848 posArray[2][ichannel] = posL[0];
1849 posArray[3][ichannel] = posL[1];
1850 posArray[4][ichannel] = posG[0];
1851 posArray[5][ichannel] = posG[1];
1852 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
1853 posArray[7][ichannel] = ichannel;
1855 // loop over array containing AliTPCCalPads
1856 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
1857 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
1858 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
1860 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
1862 (vectorArray[ivalue])[ichannel] = 0;
1867 AliTPCCalROC dummyROC(0);
1868 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
1869 AliTPCCalROC *roc = ((AliTPCCalPad*)array->At(ivalue))->GetCalROC(isector);
1870 if (!roc) roc = &dummyROC;
1871 cstream << "calPads" <<
1872 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
1873 cstream << "calPads" <<
1874 (Char_t*)((names[ivalue] + "Pad.=").Data()) << roc;
1878 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1880 cstream << "calPads" <<
1881 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
1883 cstream << "calPads" <<
1884 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
1888 cstream << "calPads" <<
1889 "sector=" << isector;
1891 cstream << "calPads" <<
1892 "row.=" << &posArray[0] <<
1893 "pad.=" << &posArray[1] <<
1894 "lx.=" << &posArray[2] <<
1895 "ly.=" << &posArray[3] <<
1896 "gx.=" << &posArray[4] <<
1897 "gy.=" << &posArray[5] <<
1898 "rpad.=" << &posArray[6] <<
1899 "channel.=" << &posArray[7];
1901 cstream << "calPads" <<
1905 delete[] vectorArray;
1906 } //for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++)
1912 delete[] mapIROCArray;
1913 delete[] mapOROCArray;
1919 void AliTPCCalibViewer::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
1921 // Write a tree with all available information
1922 // if mapFileName is speciefied, the Map information are also written to the tree
1923 // pads specified in outlierPad are not used for calculating statistics
1924 // The following statistical information on the basis of a ROC are calculated:
1925 // "_Median", "_Mean", "_LTM", "_RMS_LTM"
1926 // "_Median_OutlierCutted", "_Mean_OutlierCutted", "_RMS_OutlierCutted", "_LTM_OutlierCutted", "_RMS_LTM_OutlierCutted"
1927 // The following position variables are available:
1928 // "row", "pad", "lx", "ly", "gx", "gy", "rpad", "channel"
1930 // The tree out of this function is the basis for the AliTPCCalibViewer and the AliTPCCalibViewerGUI.
1932 AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
1934 TObjArray* mapIROCs = 0;
1935 TObjArray* mapOROCs = 0;
1936 TVectorF *mapIROCArray = 0;
1937 TVectorF *mapOROCArray = 0;
1938 Int_t mapEntries = 0;
1939 TString* mapNames = 0;
1942 TFile mapFile(mapFileName, "read");
1944 TList* listOfROCs = mapFile.GetListOfKeys();
1945 mapEntries = listOfROCs->GetEntries()/2;
1946 mapIROCs = new TObjArray(mapEntries*2);
1947 mapOROCs = new TObjArray(mapEntries*2);
1948 mapIROCArray = new TVectorF[mapEntries];
1949 mapOROCArray = new TVectorF[mapEntries];
1951 mapNames = new TString[mapEntries];
1952 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1953 TString rocName(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
1954 rocName.Remove(rocName.Length()-4, 4);
1955 mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((rocName + "IROC").Data()), ivalue);
1956 mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((rocName + "OROC").Data()), ivalue);
1957 mapNames[ivalue].Append(rocName);
1960 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1961 mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
1962 mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
1964 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
1965 (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
1966 for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
1967 (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
1970 } // if (mapFileName)
1972 TTreeSRedirector cstream(fileName);
1973 Int_t arrayEntries = 0;
1974 if (array) arrayEntries = array->GetEntries();
1976 TString* names = new TString[arrayEntries];
1977 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
1978 names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
1980 for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
1982 // get statistic for given sector
1984 TVectorF median(arrayEntries);
1985 TVectorF mean(arrayEntries);
1986 TVectorF rms(arrayEntries);
1987 TVectorF ltm(arrayEntries);
1988 TVectorF ltmrms(arrayEntries);
1989 TVectorF medianWithOut(arrayEntries);
1990 TVectorF meanWithOut(arrayEntries);
1991 TVectorF rmsWithOut(arrayEntries);
1992 TVectorF ltmWithOut(arrayEntries);
1993 TVectorF ltmrmsWithOut(arrayEntries);
1995 TVectorF *vectorArray = new TVectorF[arrayEntries];
1996 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
1997 vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
1999 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
2000 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
2001 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
2002 AliTPCCalROC* outlierROC = 0;
2003 if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
2005 median[ivalue] = calROC->GetMedian();
2006 mean[ivalue] = calROC->GetMean();
2007 rms[ivalue] = calROC->GetRMS();
2008 Double_t ltmrmsValue = 0;
2009 ltm[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction);
2010 ltmrms[ivalue] = ltmrmsValue;
2012 medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
2013 meanWithOut[ivalue] = calROC->GetMean(outlierROC);
2014 rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
2016 ltmWithOut[ivalue] = calROC->GetLTM(<mrmsValue, ltmFraction, outlierROC);
2017 ltmrmsWithOut[ivalue] = ltmrmsValue;
2021 median[ivalue] = 0.;
2025 ltmrms[ivalue] = 0.;
2026 medianWithOut[ivalue] = 0.;
2027 meanWithOut[ivalue] = 0.;
2028 rmsWithOut[ivalue] = 0.;
2029 ltmWithOut[ivalue] = 0.;
2030 ltmrmsWithOut[ivalue] = 0.;
2035 // fill vectors of variable per pad
2037 TVectorF *posArray = new TVectorF[8];
2038 for (Int_t ivalue = 0; ivalue < 8; ivalue++)
2039 posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
2041 Float_t posG[3] = {0};
2042 Float_t posL[3] = {0};
2044 for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
2045 for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
2046 tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
2047 tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
2048 posArray[0][ichannel] = irow;
2049 posArray[1][ichannel] = ipad;
2050 posArray[2][ichannel] = posL[0];
2051 posArray[3][ichannel] = posL[1];
2052 posArray[4][ichannel] = posG[0];
2053 posArray[5][ichannel] = posG[1];
2054 posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
2055 posArray[7][ichannel] = ichannel;
2057 // loop over array containing AliTPCCalPads
2058 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
2059 AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
2060 AliTPCCalROC* calROC = calPad->GetCalROC(isector);
2062 (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
2064 (vectorArray[ivalue])[ichannel] = 0;
2070 cstream << "calPads" <<
2071 "sector=" << isector;
2073 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
2074 cstream << "calPads" <<
2075 (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
2076 (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
2077 (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
2078 (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
2079 (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
2081 cstream << "calPads" <<
2082 (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
2083 (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
2084 (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
2085 (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
2086 (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
2090 for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
2091 cstream << "calPads" <<
2092 (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
2096 for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
2098 cstream << "calPads" <<
2099 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
2101 cstream << "calPads" <<
2102 (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
2106 cstream << "calPads" <<
2107 "row.=" << &posArray[0] <<
2108 "pad.=" << &posArray[1] <<
2109 "lx.=" << &posArray[2] <<
2110 "ly.=" << &posArray[3] <<
2111 "gx.=" << &posArray[4] <<
2112 "gy.=" << &posArray[5] <<
2113 "rpad.=" << &posArray[6] <<
2114 "channel.=" << &posArray[7];
2116 cstream << "calPads" <<
2120 delete[] vectorArray;
2128 delete[] mapIROCArray;
2129 delete[] mapOROCArray;
2135 void AliTPCCalibViewer::MakeTree(const char *outPutFileName, const Char_t *inputFileName, AliTPCCalPad *outlierPad, Float_t ltmFraction, const char *mapFileName ){
2137 // Function to create a calibration Tree with all available information.
2138 // See also documentation to MakeTree
2139 // the file "inputFileName" must be in the following format (see also CreateObjectList):
2140 // (each colum separated by tabs, "dependingOnType" can have 2 or 3 colums)
2142 // type path dependingOnType
2145 // dependingOnType = CETmean CEQmean CETrms
2147 // type == "Pulser":
2148 // dependingOnType = PulserTmean PulsterQmean PulserTrms
2150 // type == "Pedestals":
2151 // dependingOnType = Pedestals Noise
2153 // type == "CalPad":
2154 // dependingOnType = NameInFile NameToWriteToFile
2158 CreateObjectList(inputFileName, &objArray);
2159 MakeTree(outPutFileName, &objArray, mapFileName, outlierPad, ltmFraction);
2163 void AliTPCCalibViewer::CreateObjectList(const Char_t *filename, TObjArray *calibObjects){
2165 // Function to create a TObjArray out of a given file
2166 // the file must be in the following format:
2167 // (each colum separated by tabs, "dependingOnType" can have 2 or 3 colums)
2170 // type path dependingOnType
2173 // dependingOnType = CETmean CEQmean CETrms
2175 // type == "Pulser":
2176 // dependingOnType = PulserTmean PulsterQmean PulserTrms
2178 // type == "Pedestals":
2179 // dependingOnType = Pedestals Noise
2181 // type == "CalPad":
2182 // dependingOnType = NameInFile NameToWriteToFile
2186 if ( calibObjects == 0x0 ) return;
2189 if ( !in.is_open() ){
2190 fprintf(stderr,"Error: cannot open list file '%s'", filename);
2194 AliTPCCalPad *calPad=0x0;
2200 TObjArray *arrFileLine = sFile.Tokenize("\n");
2201 TIter nextLine(arrFileLine);
2203 TObjString *sObjLine = 0x0;
2204 while ( (sObjLine = (TObjString*)nextLine()) ){
2205 TString sLine(sObjLine->GetString());
2207 TObjArray *arrCol = sLine.Tokenize("\t");
2208 Int_t nCols = arrCol->GetEntriesFast();
2210 TObjString *sObjType = (TObjString*)(arrCol->At(0));
2211 TObjString *sObjFileName = (TObjString*)(arrCol->At(1));
2212 TObjString *sObjName = 0x0;
2214 if ( !sObjType || !sObjFileName ) continue;
2215 TString sType(sObjType->GetString());
2216 TString sFileName(sObjFileName->GetString());
2217 printf("Type %s, opening %s \n", sType.Data(), sFileName.Data());
2218 TFile *fIn = TFile::Open(sFileName);
2220 fprintf(stderr,"File not found: '%s'", sFileName.Data());
2224 if ( sType == "CE" ){ // next three colums are the names for CETmean, CEQmean and CETrms
2225 AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
2226 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());
2227 if (nCols > 2) { // check, if the name is provided
2228 sObjName = (TObjString*)(arrCol->At(2));
2229 calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2231 else calPad->SetNameTitle("CETmean","CETmean");
2232 calibObjects->Add(calPad);
2234 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());
2235 if (nCols > 3) { // check, if the name is provided
2236 sObjName = (TObjString*)(arrCol->At(3));
2237 calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2239 else calPad->SetNameTitle("CEQmean","CEQmean");
2240 calibObjects->Add(calPad);
2242 calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
2243 if (nCols > 4) { // check, if the name is provided
2244 sObjName = (TObjString*)(arrCol->At(4));
2245 calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2247 else calPad->SetNameTitle("CETrms","CETrms");
2248 calibObjects->Add(calPad);
2250 } else if ( sType == "Pulser") {
2251 AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
2253 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());
2255 sObjName = (TObjString*)(arrCol->At(2));
2256 calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2258 else calPad->SetNameTitle("PulserTmean","PulserTmean");
2259 calibObjects->Add(calPad);
2261 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());
2263 sObjName = (TObjString*)(arrCol->At(3));
2264 calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2266 else calPad->SetNameTitle("PulserQmean","PulserQmean");
2267 calibObjects->Add(calPad);
2269 calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
2271 sObjName = (TObjString*)(arrCol->At(4));
2272 calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2274 else calPad->SetNameTitle("PulserTrms","PulserTrms");
2275 calibObjects->Add(calPad);
2277 } else if ( sType == "Pedestals") {
2278 AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
2280 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());
2282 sObjName = (TObjString*)(arrCol->At(2));
2283 calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2285 else calPad->SetNameTitle("Pedestals","Pedestals");
2286 calibObjects->Add(calPad);
2288 calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());
2290 sObjName = (TObjString*)(arrCol->At(3));
2291 calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2293 else calPad->SetNameTitle("Noise","Noise");
2294 calibObjects->Add(calPad);
2296 } else if ( sType == "CalPad") {
2297 if (nCols > 2) sObjName = (TObjString*)(arrCol->At(2));
2299 calPad = new AliTPCCalPad(*(AliTPCCalPad*)fIn->Get(sObjName->GetString().Data()));
2301 sObjName = (TObjString*)(arrCol->At(3));
2302 calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2304 calibObjects->Add(calPad);
2306 fprintf(stderr,"Undefined Type: '%s'",sType.Data());