]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCalibViewer.cxx
13e5dfd4a6a92ae420e0944163eeec461d371f14
[u/mrichter/AliRoot.git] / TPC / AliTPCCalibViewer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 //  Class for viewing/visualizing TPC calibration data                       //
20 //  base on  TTree functionality for visualization                           //
21 //                                                                           //
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")    //
26 //  Have fun!                                                                //
27 //  EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")      //
28 //                                                                           //
29 //  If you like to click, we recommand you the                               //
30 //    AliTPCCalibViewerGUI                                                   //
31 //                                                                           //
32 //    THE DOCUMENTATION IS STILL NOT COMPLETED !!!!                          //
33 //                                                                           //
34 ///////////////////////////////////////////////////////////////////////////////
35
36 //
37 // ROOT includes 
38 //
39
40 #include <fstream>
41 #include <iostream>
42
43 #include <TFile.h>
44 #include <TFriendElement.h>
45 #include <TGraph.h>
46 #include <TKey.h>
47 #include <TPad.h>
48 //#include <TCanvas.h>
49 #include <TH1.h> 
50 #include <TH1F.h>
51 #include <TLegend.h>
52 #include <TLine.h>
53 #include <TMath.h>
54 #include <TObjString.h>
55 //#include <TROOT.h>
56 #include <TRandom.h>
57 #include <TString.h>
58 #include <TStyle.h>
59 #include <TTreeStream.h>
60
61 #include "AliTPCCalibCE.h"
62 #include "AliMathBase.h"
63 #include "AliTPCCalPad.h"
64 #include "AliTPCCalROC.h"
65 #include "AliTPCCalibPedestal.h"
66 #include "AliTPCCalibPulser.h"
67
68 //
69 // AliRoot includes
70 //
71 #include "AliTPCCalibViewer.h"
72
73 ClassImp(AliTPCCalibViewer)
74
75
76 AliTPCCalibViewer::AliTPCCalibViewer()
77                   :TObject(),
78                    fTree(0),
79                    fFile(0),
80                    fListOfObjectsToBeDeleted(0),
81                    fTreeMustBeDeleted(0), 
82                    fAbbreviation(0), 
83                    fAppendString(0)
84 {
85   //
86   // Default constructor
87   //
88
89 }
90
91 //_____________________________________________________________________________
92 AliTPCCalibViewer::AliTPCCalibViewer(const AliTPCCalibViewer &c)
93                   :TObject(c),
94                    fTree(0),
95                    fFile(0),
96                    fListOfObjectsToBeDeleted(0),
97                    fTreeMustBeDeleted(0),
98                    fAbbreviation(0), 
99                    fAppendString(0)
100 {
101   //
102   // dummy AliTPCCalibViewer copy constructor
103   // not yet working!!!
104   //
105   fTree = c.fTree;
106   fTreeMustBeDeleted = c.fTreeMustBeDeleted;
107   //fFile = new TFile(*(c.fFile));
108   fListOfObjectsToBeDeleted = c.fListOfObjectsToBeDeleted;
109   fAbbreviation = c.fAbbreviation;
110   fAppendString = c.fAppendString;
111 }
112
113 //_____________________________________________________________________________
114 AliTPCCalibViewer::AliTPCCalibViewer(TTree *const tree)
115                   :TObject(),
116                    fTree(0),
117                    fFile(0),
118                    fListOfObjectsToBeDeleted(0),
119                    fTreeMustBeDeleted(0),
120                    fAbbreviation(0), 
121                    fAppendString(0)
122 {
123   //
124   // Constructor that initializes the calibration viewer
125   //
126   fTree = tree;
127   fTreeMustBeDeleted = kFALSE;
128   fListOfObjectsToBeDeleted = new TObjArray();
129   fAbbreviation = "~";
130   fAppendString = ".fElements";
131 }
132
133 //_____________________________________________________________________________
134 AliTPCCalibViewer::AliTPCCalibViewer(const char* fileName, const char* treeName)
135                   :TObject(),
136                    fTree(0),
137                    fFile(0),
138                    fListOfObjectsToBeDeleted(0),
139                    fTreeMustBeDeleted(0),
140                    fAbbreviation(0), 
141                    fAppendString(0)
142                    
143 {
144    //
145    // Constructor to initialize the calibration viewer
146    // the file 'fileName' contains the tree 'treeName'
147    //
148    fFile = new TFile(fileName, "read");
149    fTree = (TTree*) fFile->Get(treeName);
150    fTreeMustBeDeleted = kTRUE;
151    fListOfObjectsToBeDeleted = new TObjArray();
152    fAbbreviation = "~";
153    fAppendString = ".fElements";
154 }
155                    
156 //____________________________________________________________________________
157 AliTPCCalibViewer & AliTPCCalibViewer::operator =(const AliTPCCalibViewer & param)
158 {
159    //
160    // assignment operator - dummy
161    // not yet working!!!
162    //
163    fTree = param.fTree;
164    fTreeMustBeDeleted = param.fTreeMustBeDeleted;
165    //fFile = new TFile(*(param.fFile));
166    fListOfObjectsToBeDeleted = param.fListOfObjectsToBeDeleted;
167    fAbbreviation = param.fAbbreviation;
168    fAppendString = param.fAppendString;
169    return (*this);
170 }
171
172 //_____________________________________________________________________________
173 AliTPCCalibViewer::~AliTPCCalibViewer()
174 {
175    //
176    // AliTPCCalibViewer destructor
177    // all objects will be deleted, the file will be closed, the pictures will disappear
178    //
179    if (fTree && fTreeMustBeDeleted) {
180       fTree->SetCacheSize(0);
181       fTree->Delete();
182       //delete fTree;
183    }
184    if (fFile) {
185       fFile->Close();
186       fFile = 0;
187    }
188
189    for (Int_t i = fListOfObjectsToBeDeleted->GetEntriesFast()-1; i >= 0; i--) {
190       //cout << "Index " << i << " trying to delete the following object: " << fListOfObjectsToBeDeleted->At(i)->GetName() << "..."<< endl;
191       delete fListOfObjectsToBeDeleted->At(i);
192    }
193    delete fListOfObjectsToBeDeleted;
194 }
195
196 //_____________________________________________________________________________
197 void AliTPCCalibViewer::Delete(Option_t* option) {
198    //
199    // Should be called from AliTPCCalibViewerGUI class only.
200    // If you use Delete() do not call the destructor.
201    // All objects (except those contained in fListOfObjectsToBeDeleted) will be deleted, the file will be closed.
202    //
203    
204    option = option;  // to avoid warnings on compiling   
205    if (fTree && fTreeMustBeDeleted) {
206       fTree->SetCacheSize(0);
207       fTree->Delete();
208    }
209    if (fFile)
210       delete fFile;
211    delete fListOfObjectsToBeDeleted;
212 }
213
214
215 const char* AliTPCCalibViewer::AddAbbreviations(const Char_t *c, Bool_t printDrawCommand){ 
216    // Replace all "<variable>" with "<variable><fAbbreviation>" (Adds forgotten "~")
217    // but take care on the statistical information, like "CEQmean_Mean"
218    // and also take care on correct given variables, like "CEQmean~"
219    // 
220    // For each variable out of "listOfVariables":
221    // - 'Save' correct items:
222    //   - form <replaceString>, take <variable>'s first char, add <removeString>, add rest of <variable>, e.g. "C!#EQmean" (<removeString> = "!#")
223    //   - For each statistical information in "listOfNormalizationVariables":
224    //     - ReplaceAll <variable><statistical_Information> with <replaceString><statistical_Information>
225    //   - ReplaceAll <variable><abbreviation> with <replaceString><abbreviation>, e.g. "CEQmean~" -> "C!#EQmean~"
226    //   - ReplaceAll <variable><appendStr> with <replaceString><appendStr>, e.g. "CEQmean.fElements" -> "C!#EQmean.fElements"
227    //
228    // - Do actual replacing:
229    //   - ReplaceAll <variable> with <variable><fAbbreviation>, e.g. "CEQmean" -> "CEQmean~"
230    //
231    // - Undo saving:
232    //   - For each statistical information in "listOfNormalizationVariables":
233    //     - ReplaceAll <replaceString><statistical_Information> with <variable><statistical_Information> 
234    //   - ReplaceAll <replaceString><abbreviation> with <variable><abbreviation>, e.g. "C!#EQmean~" -> "CEQmean~"
235    //   - ReplaceAll <replaceString><appendStr> with <variable><appendStr>, e.g. "C!#EQmean.fElements" -> "CEQmean.fElements"
236    // 
237    // Now all the missing "~" should be added.
238    
239    TString str(c);
240    TString removeString = "!#";  // very unpropable combination of chars
241    TString replaceString = "";
242    TString searchString = "";
243    TString normString = "";
244    TObjArray *listOfVariables = GetListOfVariables();
245    listOfVariables->Add(new TObjString("channel"));
246    listOfVariables->Add(new TObjString("gx"));
247    listOfVariables->Add(new TObjString("gy"));
248    listOfVariables->Add(new TObjString("lx"));
249    listOfVariables->Add(new TObjString("ly"));
250    listOfVariables->Add(new TObjString("pad"));
251    listOfVariables->Add(new TObjString("row"));
252    listOfVariables->Add(new TObjString("rpad"));
253    listOfVariables->Add(new TObjString("sector"));
254    TObjArray *listOfNormalizationVariables = GetListOfNormalizationVariables();
255    Int_t nVariables = listOfVariables->GetEntriesFast();
256    Int_t nNorm = listOfNormalizationVariables->GetEntriesFast();
257    
258    Int_t *varLengths = new Int_t[nVariables];
259    for (Int_t i = 0; i < nVariables; i++) {
260       varLengths[i] = ((TObjString*)listOfVariables->At(i))->String().Length();
261    }
262    Int_t *normLengths = new Int_t[nNorm];
263    for (Int_t i = 0; i < nNorm; i++) {
264       normLengths[i] = ((TObjString*)listOfNormalizationVariables->At(i))->String().Length();
265       // printf("normLengths[%i] (%s) = %i \n", i,((TObjString*)listOfNormalizationVariables->At(i))->String().Data(), normLengths[i]);
266    }
267    Int_t *varSort = new Int_t[nVariables];
268    TMath::Sort(nVariables, varLengths, varSort, kTRUE);
269    Int_t *normSort = new Int_t[nNorm];
270    TMath::Sort(nNorm, normLengths, normSort, kTRUE);
271    // for (Int_t i = 0; i<nNorm; i++)  printf("normLengths: %i\n", normLengths[normSort[i]]);
272    // for (Int_t i = 0; i<nVariables; i++) printf("varLengths: %i\n", varLengths[varSort[i]]);
273    
274    for (Int_t ivar = 0; ivar < nVariables; ivar++) {
275       // ***** save correct tokens *****
276       // first get the next variable:
277       searchString = ((TObjString*)listOfVariables->At(varSort[ivar]))->String();
278       // printf("searchString: %s ++++++++++++++\n", searchString.Data());
279       // form replaceString:
280       replaceString = "";
281       for (Int_t i = 0; i < searchString.Length(); i++) {
282          replaceString.Append(searchString[i]);
283          if (i == 0) replaceString.Append(removeString);
284       }
285       // go through normalization:
286       // printf("go through normalization\n");
287       for (Int_t inorm = 0; inorm < nNorm; inorm++) {
288          // printf(" inorm=%i, nNorm=%i, normSort[inorm]=%i \n", inorm, nNorm, normSort[inorm]);
289          normString = ((TObjString*)listOfNormalizationVariables->At(normSort[inorm]))->String();
290          // printf(" walking in normalization, i=%i, normString=%s \n", inorm, normString.Data());
291          str.ReplaceAll(searchString + normString, replaceString + normString);
292          // like: str.ReplaceAll("CEQmean_Mean", "C!EQmean_Mean");
293       }
294       str.ReplaceAll(searchString + fAbbreviation, replaceString + fAbbreviation);
295       // like: str.ReplaceAll("CEQmean~", "C!EQmean~");
296       str.ReplaceAll(searchString + fAppendString,    replaceString + fAppendString);
297       // like: str.ReplaceAll("CEQmean.fElements", "C!EQmean.fElements");
298       
299       // ***** add missing extensions *****
300       str.ReplaceAll(searchString, replaceString + fAbbreviation);
301       // like: str.ReplaceAll("CEQmean", "C!EQmean~");
302    }
303    
304    // ***** undo saving *****
305    str.ReplaceAll(removeString, "");
306   
307    if (printDrawCommand) std::cout << "The string looks now like: " << str.Data() << std::endl;
308    delete varSort;
309    delete normSort;
310    return str.Data();
311 }
312
313
314
315
316 //_____________________________________________________________________________
317 Int_t AliTPCCalibViewer::EasyDraw(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
318   //
319   // easy drawing of data, use '~' for abbreviation of '.fElements'
320   // example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
321  // sector: sector-number - only the specified sector will be drwawn
322   //         'A'/'C' or 'a'/'c' - side A/C will be drawn
323   //         'ALL' - whole TPC will be drawn, projected on one side
324   // cuts: specifies cuts
325   // drawOptions: draw options like 'same'
326   // writeDrawCommand: write the command, that is passed to TTree::Draw
327   //
328
329    TString drawStr(drawCommand);
330    TString sectorStr(sector);
331    sectorStr.ToUpper();
332    TString cutStr("");
333    //TString drawOptionsStr("profcolz ");
334    Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
335    if (dangerousToDraw) {
336       Warning("EasyDraw", "The draw string must not contain ':' or '>>'. Using only first variable for drawing!");
337 //      return -1;
338 //      drawStr.Resize(drawStr.First(">"));
339       drawStr.Resize(drawStr.First(":"));
340    }
341
342    TString drawOptionsStr("");
343    TRandom rnd(0);
344    Int_t rndNumber = rnd.Integer(10000);
345
346    if (drawOptions && strcmp(drawOptions, "") != 0)
347       drawOptionsStr += drawOptions;
348    else
349       drawOptionsStr += "profcolz";
350
351    if (sectorStr == "A") {
352       drawStr += Form(":gy%s:gx%s>>prof", fAppendString.Data(), fAppendString.Data());
353       drawStr += rndNumber;
354       drawStr += "(330,-250,250,330,-250,250)";
355       cutStr += "(sector/18)%2==0 ";
356    }
357    else if  (sectorStr == "C") {
358       drawStr += Form(":gy%s:gx%s>>prof", fAppendString.Data(), fAppendString.Data());
359       drawStr += rndNumber;
360       drawStr += "(330,-250,250,330,-250,250)";
361       cutStr += "(sector/18)%2==1 ";
362    }
363    else if  (sectorStr == "ALL") {
364       drawStr += Form(":gy%s:gx%s>>prof", fAppendString.Data(), fAppendString.Data());
365       drawStr += rndNumber;
366       drawStr += "(330,-250,250,330,-250,250)";
367    }
368    else if  (sectorStr.Contains("S")) {
369       drawStr += Form(":rpad%s:row%s+(sector>35)*63>>prof", fAppendString.Data(), fAppendString.Data());
370       drawStr += rndNumber;
371       drawStr += "(159,0,159,140,-70,70)";
372       TString sec=sectorStr;
373       sec.Remove(0,1);
374       cutStr += "sector%36=="+sec+" ";
375    }
376    else if (sectorStr.IsDigit()) {
377       Int_t isec = sectorStr.Atoi();
378       drawStr += Form(":rpad%s:row%s>>prof", fAppendString.Data(), fAppendString.Data());
379       drawStr += rndNumber;
380       if (isec < 36 && isec >= 0)
381          drawStr += "(63,0,63,108,-54,54)";
382       else if (isec < 72 && isec >= 36)
383          drawStr += "(96,0,96,140,-70,70)";
384       else {
385          Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
386          return -1;
387       }
388       cutStr += "(sector==";
389       cutStr += isec;
390       cutStr += ") ";
391    }
392
393    if (cuts && cuts[0] != 0) {
394       if (cutStr.Length() != 0) cutStr += "&& ";
395       cutStr += "(";
396       cutStr += cuts;
397       cutStr += ")";
398    }
399    drawStr.ReplaceAll(fAbbreviation, fAppendString);
400    cutStr.ReplaceAll(fAbbreviation, fAppendString);
401    if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" <<  cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
402    Int_t returnValue = fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
403    TString profName("prof");
404    profName += rndNumber;
405    TObject *obj = gDirectory->Get(profName.Data());
406    if (obj && obj->InheritsFrom("TH1")) FormatHistoLabels((TH1*)obj);
407    return returnValue;
408 }
409
410
411 Int_t AliTPCCalibViewer::EasyDraw(const char* drawCommand, Int_t sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
412   //
413   // easy drawing of data, use '~' for abbreviation of '.fElements'
414   // example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
415   // sector: sector-number - only the specified sector will be drwawn
416   // cuts: specifies cuts
417   // drawOptions: draw options like 'same'
418   // writeDrawCommand: write the command, that is passed to TTree::Draw
419   //
420    if (sector >= 0 && sector < 72) {
421       char sectorChr[3];
422       sprintf(sectorChr, "%i", sector);
423       return EasyDraw(drawCommand, sectorChr, cuts, drawOptions, writeDrawCommand);
424    }
425    Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
426    return -1;
427 }
428
429
430 //_____________________________________________________________________________
431 Int_t AliTPCCalibViewer::EasyDraw1D(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
432   //
433   // easy drawing of data, use '~' for abbreviation of '.fElements'
434   // example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
435   // sector: sector-number - the specified sector will be drwawn
436   //         'A'/'C' or 'a'/'c' - side A/C will be drawn
437   //         'ALL' - whole TPC will be drawn, projected on one side
438   // cuts: specifies cuts
439   // drawOptions: draw options like 'same'
440   // writeDrawCommand: write the command, that is passed to TTree::Draw
441   //
442
443    TString drawStr(drawCommand);
444    TString sectorStr(sector);
445    TString drawOptionsStr(drawOptions);
446    sectorStr.ToUpper();
447    TString cutStr("");
448
449    if (sectorStr == "A")
450       cutStr += "(sector/18)%2==0 ";
451    else if  (sectorStr == "C")
452       cutStr += "(sector/18)%2==1 ";
453    else if (sectorStr.IsDigit()) {
454       Int_t isec = sectorStr.Atoi();
455       if (isec < 0 || isec > 71) {
456          Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
457          return -1;
458       }
459       cutStr += "(sector==";
460       cutStr += isec;
461       cutStr += ") ";
462    }
463    else if  (sectorStr.Contains("S")) {
464       TString sec=sectorStr;
465       sec.Remove(0,1);
466       cutStr += "sector%36=="+sec+" ";
467    }
468
469    if (cuts && cuts[0] != 0) {
470       if (cutStr.Length() != 0) cutStr += "&& ";
471       cutStr += "(";
472       cutStr += cuts;
473       cutStr += ")";
474    }
475
476    drawStr.ReplaceAll(fAbbreviation, fAppendString);
477    cutStr.ReplaceAll(fAbbreviation, fAppendString);
478    if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" <<  cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
479    Int_t returnValue = fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
480    if (returnValue == -1) return -1;
481    
482    TObject *obj = (gPad) ? gPad->GetPrimitive("htemp") : 0; 
483    if (!obj) obj = (TH1F*)gDirectory->Get("htemp");
484    if (!obj) obj = gPad->GetPrimitive("tempHist");
485    if (!obj) obj = (TH1F*)gDirectory->Get("tempHist");
486    if (!obj) obj = gPad->GetPrimitive("Graph");
487    if (!obj) obj = (TH1F*)gDirectory->Get("Graph");
488    if (obj && obj->InheritsFrom("TH1")) FormatHistoLabels((TH1*)obj);
489    return returnValue;
490 }
491
492
493 Int_t AliTPCCalibViewer::EasyDraw1D(const char* drawCommand, Int_t sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
494   //
495   // easy drawing of data, use '~' for abbreviation of '.fElements'
496   // example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
497   // sector: sector-number - the specified sector will be drwawn
498   // cuts: specifies cuts
499   // drawOptions: draw options like 'same'
500   // writeDrawCommand: write the command, that is passed to TTree::Draw
501   //
502
503    if (sector >= 0 && sector < 72) {
504       char sectorChr[3];
505       sprintf(sectorChr, "%i", sector);
506       return EasyDraw1D(drawCommand, sectorChr, cuts, drawOptions, writeDrawCommand);
507    }
508   Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
509   return -1;
510 }
511
512
513 void AliTPCCalibViewer::FormatHistoLabels(TH1 *histo) const {
514    // 
515    // formats title and axis labels of histo 
516    // removes '.fElements'
517    // 
518    if (!histo) return;
519    TString replaceString(fAppendString.Data());
520    TString *str = new TString(histo->GetTitle());
521    str->ReplaceAll(replaceString, "");
522    histo->SetTitle(str->Data());
523    delete str;
524    if (histo->GetXaxis()) {
525       str = new TString(histo->GetXaxis()->GetTitle());
526       str->ReplaceAll(replaceString, "");
527       histo->GetXaxis()->SetTitle(str->Data());
528       delete str;
529    }
530    if (histo->GetYaxis()) {
531       str = new TString(histo->GetYaxis()->GetTitle());
532       str->ReplaceAll(replaceString, "");
533       histo->GetYaxis()->SetTitle(str->Data());
534       delete str;
535    }
536    if (histo->GetZaxis()) {
537       str = new TString(histo->GetZaxis()->GetTitle());
538       str->ReplaceAll(replaceString, "");
539       histo->GetZaxis()->SetTitle(str->Data());
540       delete str;
541    }
542 }
543
544
545 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 {
546    // 
547    // Easy drawing of data, in principle the same as EasyDraw1D
548    // Difference: A line for the mean / median / LTM is drawn 
549    // in 'sigmas' you can specify in which distance to the mean/median/LTM you want to see a line in sigma-units, separated by ';'
550    // 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.
551    // "plotMean", "plotMedian" and "plotLTM": what kind of lines do you want to see?
552    // 
553    if (sector >= 0 && sector < 72) {
554       char sectorChr[3];
555       sprintf(sectorChr, "%i", sector);
556       return DrawHisto1D(drawCommand, sectorChr, cuts, sigmas, plotMean, plotMedian, plotLTM);
557    }
558    Error("DrawHisto1D","The TPC contains only sectors between 0 and 71.");
559    return -1;
560 }   
561
562
563 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 {
564    // 
565    // Easy drawing of data, in principle the same as EasyDraw1D
566    // Difference: A line for the mean / median / LTM is drawn 
567    // in 'sigmas' you can specify in which distance to the mean/median/LTM you want to see a line in sigma-units, separated by ';'
568    // 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.
569    // "plotMean", "plotMedian" and "plotLTM": what kind of lines do you want to see?
570    // 
571    Int_t oldOptStat = gStyle->GetOptStat();
572    gStyle->SetOptStat(0000000);
573    Double_t ltmFraction = 0.8;
574    
575    TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");  
576    TVectorF nsigma(sigmasTokens->GetEntriesFast());
577    for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
578       TString str(((TObjString*)sigmasTokens->At(i))->GetString());
579       Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
580       nsigma[i] = sig;
581    }
582    
583    TString drawStr(drawCommand);
584    Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
585    if (dangerousToDraw) {
586       Warning("DrawHisto1D", "The draw string must not contain ':' or '>>'.");
587       return -1;
588    }
589    drawStr += " >> tempHist";
590    Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts);
591    TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
592    // FIXME is this histogram deleted automatically?
593    Double_t *values = fTree->GetV1();  // value is the array containing 'entries' numbers
594    
595    Double_t mean = TMath::Mean(entries, values);
596    Double_t median = TMath::Median(entries, values);
597    Double_t sigma = TMath::RMS(entries, values);
598    Double_t maxY = htemp->GetMaximum();
599    
600    char c[500];
601    TLegend * legend = new TLegend(.7,.7, .99, .99, "Statistical information");
602 //    sprintf(c, "%s, sector: %i", type, sector);
603    //fListOfObjectsToBeDeleted->Add(legend);
604
605    if (plotMean) {
606       // draw Mean
607       TLine* line = new TLine(mean, 0, mean, maxY);
608       //fListOfObjectsToBeDeleted->Add(line);
609       line->SetLineColor(kRed);
610       line->SetLineWidth(2);
611       line->SetLineStyle(1);
612       line->Draw();
613       sprintf(c, "Mean: %f", mean);
614       legend->AddEntry(line, c, "l");
615       // draw sigma lines
616       for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
617          TLine* linePlusSigma = new TLine(mean + nsigma[i] * sigma, 0, mean + nsigma[i] * sigma, maxY);
618          //fListOfObjectsToBeDeleted->Add(linePlusSigma);
619          linePlusSigma->SetLineColor(kRed);
620          linePlusSigma->SetLineStyle(2 + i);
621          linePlusSigma->Draw();
622          TLine* lineMinusSigma = new TLine(mean - nsigma[i] * sigma, 0, mean - nsigma[i] * sigma, maxY);
623          //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
624          lineMinusSigma->SetLineColor(kRed);
625          lineMinusSigma->SetLineStyle(2 + i);
626          lineMinusSigma->Draw();
627          sprintf(c, "%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i] * sigma));
628          legend->AddEntry(lineMinusSigma, c, "l");
629       }
630    }
631    if (plotMedian) {
632       // draw median
633       TLine* line = new TLine(median, 0, median, maxY);
634       //fListOfObjectsToBeDeleted->Add(line);
635       line->SetLineColor(kBlue);
636       line->SetLineWidth(2);
637       line->SetLineStyle(1);
638       line->Draw();
639       sprintf(c, "Median: %f", median);
640       legend->AddEntry(line, c, "l");
641       // draw sigma lines
642       for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
643          TLine* linePlusSigma = new TLine(median + nsigma[i] * sigma, 0, median + nsigma[i]*sigma, maxY);
644          //fListOfObjectsToBeDeleted->Add(linePlusSigma);
645          linePlusSigma->SetLineColor(kBlue);
646          linePlusSigma->SetLineStyle(2 + i);
647          linePlusSigma->Draw();
648          TLine* lineMinusSigma = new TLine(median - nsigma[i] * sigma, 0, median - nsigma[i]*sigma, maxY);
649          //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
650          lineMinusSigma->SetLineColor(kBlue);
651          lineMinusSigma->SetLineStyle(2 + i);
652          lineMinusSigma->Draw();
653          sprintf(c, "%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i] * sigma));
654          legend->AddEntry(lineMinusSigma, c, "l");
655       }
656    }
657    if (plotLTM) {
658       // draw LTM
659       Double_t ltmRms = 0;
660       Double_t ltm = GetLTM(entries, values, &ltmRms, ltmFraction);
661       TLine* line = new TLine(ltm, 0, ltm, maxY);
662       //fListOfObjectsToBeDeleted->Add(line);
663       line->SetLineColor(kGreen+2);
664       line->SetLineWidth(2);
665       line->SetLineStyle(1);
666       line->Draw();
667       sprintf(c, "LTM: %f", ltm);
668       legend->AddEntry(line, c, "l");
669       // draw sigma lines
670       for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
671          TLine* linePlusSigma = new TLine(ltm + nsigma[i] * ltmRms, 0, ltm + nsigma[i] * ltmRms, maxY);
672          //fListOfObjectsToBeDeleted->Add(linePlusSigma);
673          linePlusSigma->SetLineColor(kGreen+2);
674          linePlusSigma->SetLineStyle(2+i);
675          linePlusSigma->Draw();
676    
677          TLine* lineMinusSigma = new TLine(ltm - nsigma[i] * ltmRms, 0, ltm - nsigma[i] * ltmRms, maxY);
678          //fListOfObjectsToBeDeleted->Add(lineMinusSigma);
679          lineMinusSigma->SetLineColor(kGreen+2);
680          lineMinusSigma->SetLineStyle(2+i);
681          lineMinusSigma->Draw();
682          sprintf(c, "%i #sigma = %f", (Int_t)(nsigma[i]), (Float_t)(nsigma[i] * ltmRms));
683          legend->AddEntry(lineMinusSigma, c, "l");
684       }
685    }
686    if (!plotMean && !plotMedian && !plotLTM) return -1;
687    legend->Draw();
688    gStyle->SetOptStat(oldOptStat);
689    return 1;
690 }
691
692
693 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 {
694    //
695    // 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
696    // The data of the distribution Begin_Latex f(x, #mu, #sigma) End_Latex are given in 'array', 'n' specifies the length of the array
697    // '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
698    // 'nbins': number of bins, 'binLow': first bin, 'binUp': last bin
699    // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, Begin_Latex t #sigma End_Latex)
700    // sigmaStep: the binsize of the generated histogram
701    // Begin_Latex 
702    // 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 }
703    // End_Latex
704    // 
705    //
706    // Creates a histogram, where you can see, how much of the data are inside sigma-intervals 
707    // around the mean/median/LTM
708    // with drawCommand, sector and cuts you specify your input data, see EasyDraw
709    // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma)
710    // sigmaStep: the binsize of the generated histogram
711    // plotMean/plotMedian/plotLTM: specifies where to put the center
712    //
713    if (sector >= 0 && sector < 72) {
714       char sectorChr[3];
715       sprintf(sectorChr, "%i", sector);
716       return SigmaCut(drawCommand, sectorChr, cuts, sigmaMax, plotMean, plotMedian, plotLTM, pm, sigmas, sigmaStep);
717    }
718    Error("SigmaCut","The TPC contains only sectors between 0 and 71.");
719    return -1;
720 }
721
722
723 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 {
724    //
725    // Creates a histogram, where you can see, how much of the data are inside sigma-intervals 
726    // around the mean/median/LTM
727    // with drawCommand, sector and cuts you specify your input data, see EasyDraw
728    // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma)
729    // sigmaStep: the binsize of the generated histogram
730    // plotMean/plotMedian/plotLTM: specifies where to put the center
731    //
732   
733    Double_t ltmFraction = 0.8;
734    
735    TString drawStr(drawCommand);
736    Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
737    if (dangerousToDraw) {
738       Warning("SigmaCut", "The draw string must not contain ':' or '>>'.");
739       return -1;
740    }
741    drawStr += " >> tempHist";
742    
743    Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
744    TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
745    // FIXME is this histogram deleted automatically?
746    Double_t *values = fTree->GetV1();  // value is the array containing 'entries' numbers
747    
748    Double_t mean = TMath::Mean(entries, values);
749    Double_t median = TMath::Median(entries, values);
750    Double_t sigma = TMath::RMS(entries, values);
751    
752    TLegend * legend = new TLegend(.7,.7, .99, .99, "Cumulative");
753    //fListOfObjectsToBeDeleted->Add(legend);
754    TH1F *cutHistoMean = 0;
755    TH1F *cutHistoMedian = 0;
756    TH1F *cutHistoLTM = 0;
757    
758    TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");  
759    TVectorF nsigma(sigmasTokens->GetEntriesFast());
760    for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
761       TString str(((TObjString*)sigmasTokens->At(i))->GetString());
762       Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
763       nsigma[i] = sig;
764    }
765   
766    if (plotMean) {
767       cutHistoMean = AliTPCCalibViewer::SigmaCut(htemp, mean, sigma, sigmaMax, sigmaStep, pm);
768       if (cutHistoMean) {
769          //fListOfObjectsToBeDeleted->Add(cutHistoMean);
770          cutHistoMean->SetLineColor(kRed);
771          legend->AddEntry(cutHistoMean, "Mean", "l");
772          cutHistoMean->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
773          cutHistoMean->Draw();
774          DrawLines(cutHistoMean, nsigma, legend, kRed, pm);
775       } // if (cutHistoMean)
776        
777    }
778    if (plotMedian) {
779       cutHistoMedian = AliTPCCalibViewer::SigmaCut(htemp, median, sigma, sigmaMax, sigmaStep, pm);
780       if (cutHistoMedian) {
781          //fListOfObjectsToBeDeleted->Add(cutHistoMedian);
782          cutHistoMedian->SetLineColor(kBlue);
783          legend->AddEntry(cutHistoMedian, "Median", "l");
784          cutHistoMedian->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
785          if (plotMean && cutHistoMean) cutHistoMedian->Draw("same");
786             else cutHistoMedian->Draw();
787          DrawLines(cutHistoMedian, nsigma, legend, kBlue, pm);
788       }  // if (cutHistoMedian)
789    }
790    if (plotLTM) {
791       Double_t ltmRms = 0;
792       Double_t ltm = GetLTM(entries, values, &ltmRms, ltmFraction);
793       cutHistoLTM = AliTPCCalibViewer::SigmaCut(htemp, ltm, ltmRms, sigmaMax, sigmaStep, pm);
794       if (cutHistoLTM) {
795          //fListOfObjectsToBeDeleted->Add(cutHistoLTM);
796          cutHistoLTM->SetLineColor(kGreen+2);
797          legend->AddEntry(cutHistoLTM, "LTM", "l");
798          cutHistoLTM->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
799          if ((plotMean && cutHistoMean) || (plotMedian && cutHistoMedian)) cutHistoLTM->Draw("same");
800             else cutHistoLTM->Draw();
801          DrawLines(cutHistoLTM, nsigma, legend, kGreen+2, pm);
802       }
803    }
804    if (!plotMean && !plotMedian && !plotLTM) return -1;
805    legend->Draw();
806    return 1;
807 }
808
809
810 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 {
811    //
812    // Creates a histogram, where you can see, how much of the data are inside sigma-intervals 
813    // around the mean/median/LTM
814    // with drawCommand, sector and cuts you specify your input data, see EasyDraw
815    // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma)
816    // sigmaStep: the binsize of the generated histogram
817    // plotMean/plotMedian/plotLTM: specifies where to put the center
818    //
819   
820    // Double_t ltmFraction = 0.8;  //unused
821    // avoid compiler warnings:
822    sigmaMax = sigmaMax;
823    pm = pm;
824    sigmaStep = sigmaStep;
825    
826    TString drawStr(drawCommand);
827    drawStr += " >> tempHist";
828    
829    Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
830    TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
831    TGraph *cutGraphMean   = 0;
832    // TGraph *cutGraphMedian = 0;
833    // TGraph *cutGraphLTM    = 0;
834    Double_t *values = fTree->GetV1();  // value is the array containing 'entries' numbers
835    Int_t    *index  = new Int_t[entries];
836    Float_t  *xarray = new Float_t[entries];
837    Float_t  *yarray = new Float_t[entries];
838    TMath::Sort(entries, values, index, kFALSE);
839    
840    Double_t mean = TMath::Mean(entries, values);
841    // Double_t median = TMath::Median(entries, values);
842    Double_t sigma = TMath::RMS(entries, values);
843    
844    TLegend * legend = new TLegend(.7,.7, .99, .99, "Cumulative");
845    //fListOfObjectsToBeDeleted->Add(legend);
846    
847    // parse sigmas string
848    TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");  
849    TVectorF nsigma(sigmasTokens->GetEntriesFast());
850    for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
851       TString str(((TObjString*)sigmasTokens->At(i))->GetString());
852       Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
853       nsigma[i] = sig;
854    }
855    
856    if (plotMean) {
857       for (Int_t i = 0; i < entries; i++) {
858          xarray[i] = TMath::Abs(values[index[i]] - mean) / sigma; 
859          yarray[i] = float(i) / float(entries);
860       }
861       cutGraphMean = new TGraph(entries, xarray, yarray);
862       if (cutGraphMean) {
863          //fListOfObjectsToBeDeleted->Add(cutGraphMean);
864          cutGraphMean->SetLineColor(kRed);
865          legend->AddEntry(cutGraphMean, "Mean", "l");
866          cutGraphMean->SetTitle(Form("%s, Cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
867          cutGraphMean->Draw("alu");
868          DrawLines(cutGraphMean, nsigma, legend, kRed, kTRUE);
869       }
870    }
871    /*
872    if (plotMedian) {
873       cutHistoMedian = AliTPCCalibViewer::SigmaCut(htemp, median, sigma, sigmaMax, sigmaStep, pm);
874       if (cutHistoMedian) {
875          fListOfObjectsToBeDeleted->Add(cutHistoMedian);
876          cutHistoMedian->SetLineColor(kBlue);
877          legend->AddEntry(cutHistoMedian, "Median", "l");
878          cutHistoMedian->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
879          if (plotMean && cutHistoMean) cutHistoMedian->Draw("same");
880             else cutHistoMedian->Draw();
881          DrawLines(cutHistoMedian, nsigma, legend, kBlue, pm);
882       }  // if (cutHistoMedian)
883    }
884    if (plotLTM) {
885       Double_t ltmRms = 0;
886       Double_t ltm = GetLTM(entries, values, &ltmRms, ltmFraction);
887       cutHistoLTM = AliTPCCalibViewer::SigmaCut(htemp, ltm, ltmRms, sigmaMax, sigmaStep, pm);
888       if (cutHistoLTM) {
889          fListOfObjectsToBeDeleted->Add(cutHistoLTM);
890          cutHistoLTM->SetLineColor(kGreen+2);
891          legend->AddEntry(cutHistoLTM, "LTM", "l");
892          cutHistoLTM->SetTitle(Form("%s, cumulative; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
893          if (plotMean && cutHistoMean || plotMedian && cutHistoMedian) cutHistoLTM->Draw("same");
894             else cutHistoLTM->Draw();
895          DrawLines(cutHistoLTM, nsigma, legend, kGreen+2, pm);
896       }
897    }*/
898    if (!plotMean && !plotMedian && !plotLTM) return -1;
899    legend->Draw();
900    return 1;
901 }
902
903
904 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 {
905    //
906    // 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"   
907    // "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
908    // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate 
909    // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
910    // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
911    // The actual work is done on the array.
912    /* Begin_Latex 
913          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 }
914       End_Latex  
915    */
916    if (sector >= 0 && sector < 72) {
917       char sectorChr[3];
918       sprintf(sectorChr, "%i", sector);
919       return Integrate(drawCommand, sectorChr, cuts, sigmaMax, plotMean, plotMedian, plotLTM, sigmas, sigmaStep);
920    }
921    Error("Integrate","The TPC contains only sectors between 0 and 71.");
922    return -1;
923    
924 }
925
926
927 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 {
928    //
929    // 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"   
930    // "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
931    // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate 
932    // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
933    // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
934    // The actual work is done on the array.
935    /* Begin_Latex 
936          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 }
937       End_Latex  
938    */
939    
940    Double_t ltmFraction = 0.8;
941    
942    TString drawStr(drawCommand);
943    drawStr += " >> tempHist";
944    
945    Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
946    TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
947    // FIXME is this histogram deleted automatically?
948    Double_t *values = fTree->GetV1();  // value is the array containing 'entries' numbers
949    
950    Double_t mean = TMath::Mean(entries, values);
951    Double_t median = TMath::Median(entries, values);
952    Double_t sigma = TMath::RMS(entries, values);
953     
954    TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");  
955    TVectorF nsigma(sigmasTokens->GetEntriesFast());
956    for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
957       TString str(((TObjString*)sigmasTokens->At(i))->GetString());
958       Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
959       nsigma[i] = sig;
960    }
961   
962    TLegend * legend = new TLegend(.7,.7, .99, .99, "Integrated histogram");
963    //fListOfObjectsToBeDeleted->Add(legend);
964    TH1F *integralHistoMean = 0;
965    TH1F *integralHistoMedian = 0;
966    TH1F *integralHistoLTM = 0;
967   
968    if (plotMean) {
969       integralHistoMean = AliTPCCalibViewer::Integrate(htemp, mean, sigma, sigmaMax, sigmaStep);
970       if (integralHistoMean) {
971          //fListOfObjectsToBeDeleted->Add(integralHistoMean);
972          integralHistoMean->SetLineColor(kRed);
973          legend->AddEntry(integralHistoMean, "Mean", "l");
974          integralHistoMean->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
975          integralHistoMean->Draw();
976          DrawLines(integralHistoMean, nsigma, legend, kRed, kTRUE);
977       }
978    }
979    if (plotMedian) {
980       integralHistoMedian = AliTPCCalibViewer::Integrate(htemp, median, sigma, sigmaMax, sigmaStep);
981       if (integralHistoMedian) {
982          //fListOfObjectsToBeDeleted->Add(integralHistoMedian);
983          integralHistoMedian->SetLineColor(kBlue);
984          legend->AddEntry(integralHistoMedian, "Median", "l");
985          integralHistoMedian->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
986          if (plotMean && integralHistoMean) integralHistoMedian->Draw("same");
987             else integralHistoMedian->Draw();
988          DrawLines(integralHistoMedian, nsigma, legend, kBlue, kTRUE);
989       }
990    }
991    if (plotLTM) {
992       Double_t ltmRms = 0;
993       Double_t ltm = GetLTM(entries, values, &ltmRms, ltmFraction);
994       integralHistoLTM = AliTPCCalibViewer::Integrate(htemp, ltm, ltmRms, sigmaMax, sigmaStep);
995       if (integralHistoLTM) {
996          //fListOfObjectsToBeDeleted->Add(integralHistoLTM);
997          integralHistoLTM->SetLineColor(kGreen+2);
998          legend->AddEntry(integralHistoLTM, "LTM", "l");
999          integralHistoLTM->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
1000          if ((plotMean && integralHistoMean) || (plotMedian && integralHistoMedian)) integralHistoLTM->Draw("same");
1001             else integralHistoLTM->Draw();
1002          DrawLines(integralHistoLTM, nsigma, legend, kGreen+2, kTRUE);
1003       }
1004    }
1005    if (!plotMean && !plotMedian && !plotLTM) return -1;
1006    legend->Draw();
1007    return 1;
1008 }
1009
1010
1011 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 {
1012    //
1013    // 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"   
1014    // "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
1015    // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate 
1016    // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
1017    // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
1018    // The actual work is done on the array.
1019    /* Begin_Latex 
1020          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 }
1021       End_Latex  
1022    */
1023    
1024    Double_t ltmFraction = 0.8;
1025    // avoid compiler warnings:
1026    sigmaMax = sigmaMax;
1027    sigmaStep = sigmaStep;
1028    
1029    TString drawStr(drawCommand);
1030    Bool_t dangerousToDraw = drawStr.Contains(":") || drawStr.Contains(">>");
1031    if (dangerousToDraw) {
1032       Warning("Integrate", "The draw string must not contain ':' or '>>'.");
1033       return -1;
1034    }
1035    drawStr += " >> tempHist";
1036    
1037    Int_t entries = EasyDraw1D(drawStr.Data(), sector, cuts, "goff");
1038    TH1F *htemp = (TH1F*)gDirectory->Get("tempHist");
1039    TGraph *integralGraphMean   = 0;
1040    TGraph *integralGraphMedian = 0;
1041    TGraph *integralGraphLTM    = 0;
1042    Double_t *values = fTree->GetV1();  // value is the array containing 'entries' numbers
1043    Int_t    *index  = new Int_t[entries];
1044    Float_t  *xarray = new Float_t[entries];
1045    Float_t  *yarray = new Float_t[entries];
1046    TMath::Sort(entries, values, index, kFALSE);
1047    
1048    Double_t mean = TMath::Mean(entries, values);
1049    Double_t median = TMath::Median(entries, values);
1050    Double_t sigma = TMath::RMS(entries, values);
1051    
1052    // parse sigmas string
1053    TObjArray *sigmasTokens = TString(sigmas).Tokenize(";");  
1054    TVectorF nsigma(sigmasTokens->GetEntriesFast());
1055    for (Int_t i = 0; i < sigmasTokens->GetEntriesFast(); i++) {
1056       TString str(((TObjString*)sigmasTokens->At(i))->GetString());
1057       Double_t sig = (str.IsFloat()) ? str.Atof() : 0;
1058       nsigma[i] = sig;
1059    }
1060   
1061    TLegend * legend = new TLegend(.7,.7, .99, .99, "Integrated histogram");
1062    //fListOfObjectsToBeDeleted->Add(legend);
1063   
1064    if (plotMean) {
1065       for (Int_t i = 0; i < entries; i++) {
1066          xarray[i] = (values[index[i]] - mean) / sigma; 
1067          yarray[i] = float(i) / float(entries);
1068       }
1069       integralGraphMean = new TGraph(entries, xarray, yarray);
1070       if (integralGraphMean) {
1071          //fListOfObjectsToBeDeleted->Add(integralGraphMean);
1072          integralGraphMean->SetLineColor(kRed);
1073          legend->AddEntry(integralGraphMean, "Mean", "l");
1074          integralGraphMean->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
1075          integralGraphMean->Draw("alu");
1076          DrawLines(integralGraphMean, nsigma, legend, kRed, kTRUE);
1077       }
1078    }
1079    if (plotMedian) {
1080       for (Int_t i = 0; i < entries; i++) {
1081          xarray[i] = (values[index[i]] - median) / sigma; 
1082          yarray[i] = float(i) / float(entries);
1083       }
1084       integralGraphMedian = new TGraph(entries, xarray, yarray);
1085       if (integralGraphMedian) {
1086          //fListOfObjectsToBeDeleted->Add(integralGraphMedian);
1087          integralGraphMedian->SetLineColor(kBlue);
1088          legend->AddEntry(integralGraphMedian, "Median", "l");
1089          integralGraphMedian->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
1090          if (plotMean && integralGraphMean) integralGraphMedian->Draw("samelu");
1091             else integralGraphMedian->Draw("alu");
1092          DrawLines(integralGraphMedian, nsigma, legend, kBlue, kTRUE);
1093       }
1094    }
1095    if (plotLTM) {
1096       Double_t ltmRms = 0;
1097       Double_t ltm = GetLTM(entries, values, &ltmRms, ltmFraction);
1098       for (Int_t i = 0; i < entries; i++) {
1099          xarray[i] = (values[index[i]] - ltm) / ltmRms; 
1100          yarray[i] = float(i) / float(entries);
1101       }
1102       integralGraphLTM = new TGraph(entries, xarray, yarray);
1103       if (integralGraphLTM) {
1104          //fListOfObjectsToBeDeleted->Add(integralGraphLTM);
1105          integralGraphLTM->SetLineColor(kGreen+2);
1106          legend->AddEntry(integralGraphLTM, "LTM", "l");
1107          integralGraphLTM->SetTitle(Form("%s, integrated; Multiples of #sigma; Fraction of included data", htemp->GetTitle()));
1108          if ((plotMean && integralGraphMean) || (plotMedian && integralGraphMedian)) integralGraphLTM->Draw("samelu");
1109             else integralGraphLTM->Draw("alu");
1110          DrawLines(integralGraphLTM, nsigma, legend, kGreen+2, kTRUE);
1111       }
1112    }
1113    if (!plotMean && !plotMedian && !plotLTM) return -1;
1114    legend->Draw();
1115    return entries;
1116 }
1117
1118
1119 void AliTPCCalibViewer::DrawLines(TH1F *histogram, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const {
1120    // 
1121    // Private function for SigmaCut(...) and Integrate(...)
1122    // Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
1123    // 
1124    
1125    // start to draw the lines, loop over requested sigmas
1126    char c[500];
1127    for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
1128       if (!pm) { 
1129          Int_t bin = histogram->GetXaxis()->FindBin(nsigma[i]);
1130          TLine* lineUp = new TLine(nsigma[i], 0, nsigma[i], histogram->GetBinContent(bin));
1131          //fListOfObjectsToBeDeleted->Add(lineUp);
1132          lineUp->SetLineColor(color);
1133          lineUp->SetLineStyle(2 + i);
1134          lineUp->Draw();
1135          TLine* lineLeft = new TLine(nsigma[i], histogram->GetBinContent(bin), 0, histogram->GetBinContent(bin));
1136          //fListOfObjectsToBeDeleted->Add(lineLeft);
1137          lineLeft->SetLineColor(color);
1138          lineLeft->SetLineStyle(2 + i);
1139          lineLeft->Draw();
1140          sprintf(c, "Fraction(%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin));
1141          legend->AddEntry(lineLeft, c, "l");
1142       }
1143       else { // if (pm)
1144          Int_t bin = histogram->GetXaxis()->FindBin(nsigma[i]);
1145          TLine* lineUp1 = new TLine(nsigma[i], 0, nsigma[i], histogram->GetBinContent(bin));
1146          //fListOfObjectsToBeDeleted->Add(lineUp1);
1147          lineUp1->SetLineColor(color);
1148          lineUp1->SetLineStyle(2 + i);
1149          lineUp1->Draw();
1150          TLine* lineLeft1 = new TLine(nsigma[i], histogram->GetBinContent(bin), histogram->GetBinLowEdge(0)+histogram->GetBinWidth(0), histogram->GetBinContent(bin));
1151          //fListOfObjectsToBeDeleted->Add(lineLeft1);
1152          lineLeft1->SetLineColor(color);
1153          lineLeft1->SetLineStyle(2 + i);
1154          lineLeft1->Draw();
1155          sprintf(c, "Fraction(+%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin));
1156          legend->AddEntry(lineLeft1, c, "l");
1157          bin = histogram->GetXaxis()->FindBin(-nsigma[i]);
1158          TLine* lineUp2 = new TLine(-nsigma[i], 0, -nsigma[i], histogram->GetBinContent(bin));
1159          //fListOfObjectsToBeDeleted->Add(lineUp2);
1160          lineUp2->SetLineColor(color);
1161          lineUp2->SetLineStyle(2 + i);
1162          lineUp2->Draw();
1163          TLine* lineLeft2 = new TLine(-nsigma[i], histogram->GetBinContent(bin), histogram->GetBinLowEdge(0)+histogram->GetBinWidth(0), histogram->GetBinContent(bin));
1164          //fListOfObjectsToBeDeleted->Add(lineLeft2);
1165          lineLeft2->SetLineColor(color);
1166          lineLeft2->SetLineStyle(2 + i);
1167          lineLeft2->Draw();
1168          sprintf(c, "Fraction(-%f #sigma) = %f",nsigma[i], histogram->GetBinContent(bin));
1169          legend->AddEntry(lineLeft2, c, "l");
1170       }
1171    }  // for (Int_t i = 0; i < nsigma.GetNoElements(); i++)   
1172 }
1173
1174
1175 void AliTPCCalibViewer::DrawLines(TGraph *graph, TVectorF nsigma, TLegend *legend, Int_t color, Bool_t pm) const {
1176    // 
1177    // Private function for SigmaCut(...) and Integrate(...)
1178    // Draws lines into the given histogram, specified by "nsigma", the lines are addeed to the legend
1179    // 
1180    
1181    // start to draw the lines, loop over requested sigmas
1182    char c[500];
1183    for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
1184       if (!pm) { 
1185          TLine* lineUp = new TLine(nsigma[i], 0, nsigma[i], graph->Eval(nsigma[i]));
1186          //fListOfObjectsToBeDeleted->Add(lineUp);
1187          lineUp->SetLineColor(color);
1188          lineUp->SetLineStyle(2 + i);
1189          lineUp->Draw();
1190          TLine* lineLeft = new TLine(nsigma[i], graph->Eval(nsigma[i]), 0, graph->Eval(nsigma[i]));
1191          //fListOfObjectsToBeDeleted->Add(lineLeft);
1192          lineLeft->SetLineColor(color);
1193          lineLeft->SetLineStyle(2 + i);
1194          lineLeft->Draw();
1195          sprintf(c, "Fraction(%f #sigma) = %f",nsigma[i], graph->Eval(nsigma[i]));
1196          legend->AddEntry(lineLeft, c, "l");
1197       }
1198       else { // if (pm)
1199          TLine* lineUp1 = new TLine(nsigma[i], 0, nsigma[i], graph->Eval(nsigma[i]));
1200          //fListOfObjectsToBeDeleted->Add(lineUp1);
1201          lineUp1->SetLineColor(color);
1202          lineUp1->SetLineStyle(2 + i);
1203          lineUp1->Draw();
1204          TLine* lineLeft1 = new TLine(nsigma[i], graph->Eval(nsigma[i]), graph->GetHistogram()->GetXaxis()->GetBinLowEdge(0), graph->Eval(nsigma[i]));
1205          //fListOfObjectsToBeDeleted->Add(lineLeft1);
1206          lineLeft1->SetLineColor(color);
1207          lineLeft1->SetLineStyle(2 + i);
1208          lineLeft1->Draw();
1209          sprintf(c, "Fraction(+%f #sigma) = %f",nsigma[i], graph->Eval(nsigma[i]));
1210          legend->AddEntry(lineLeft1, c, "l");
1211          TLine* lineUp2 = new TLine(-nsigma[i], 0, -nsigma[i], graph->Eval(-nsigma[i]));
1212          //fListOfObjectsToBeDeleted->Add(lineUp2);
1213          lineUp2->SetLineColor(color);
1214          lineUp2->SetLineStyle(2 + i);
1215          lineUp2->Draw();
1216          TLine* lineLeft2 = new TLine(-nsigma[i], graph->Eval(-nsigma[i]), graph->GetHistogram()->GetXaxis()->GetBinLowEdge(0), graph->Eval(-nsigma[i]));
1217          //fListOfObjectsToBeDeleted->Add(lineLeft2);
1218          lineLeft2->SetLineColor(color);
1219          lineLeft2->SetLineStyle(2 + i);
1220          lineLeft2->Draw();
1221          sprintf(c, "Fraction(-%f #sigma) = %f",nsigma[i], graph->Eval(-nsigma[i]));
1222          legend->AddEntry(lineLeft2, c, "l");
1223       }
1224    }  // for (Int_t i = 0; i < nsigma.GetNoElements(); i++)   
1225 }
1226
1227
1228
1229
1230
1231 /////////////////
1232 // Array tools //
1233 /////////////////
1234
1235
1236 Int_t AliTPCCalibViewer::GetBin(Float_t value, Int_t nbins, Double_t binLow, Double_t binUp){
1237    // Returns the 'bin' for 'value'
1238    // The interval between 'binLow' and 'binUp' is divided into 'nbins' equidistant bins
1239    // avoid index out of bounds error: 'if (bin < binLow) bin = binLow' and vice versa
1240    /* Begin_Latex
1241          GetBin(value) = #frac{nbins - 1}{binUp - binLow} #upoint (value - binLow) +1
1242       End_Latex
1243    */
1244    
1245    Int_t bin =  TMath::Nint( (Float_t)(value - binLow) / (Float_t)(binUp - binLow) * (nbins-1) ) + 1;
1246    // avoid index out of bounds:   
1247    if (value < binLow) bin = 0;
1248    if (value > binUp)  bin = nbins + 1;
1249    return bin;
1250    
1251 }   
1252
1253
1254 Double_t AliTPCCalibViewer::GetLTM(Int_t n, const Double_t *const array, Double_t *const sigma, Double_t fraction){
1255    //
1256    //  returns the LTM and sigma
1257    //
1258    Double_t *ddata = new Double_t[n];
1259    Double_t mean = 0, lsigma = 0;
1260    UInt_t nPoints = 0;
1261    for (UInt_t i = 0; i < (UInt_t)n; i++) {
1262          ddata[nPoints]= array[nPoints];
1263          nPoints++;
1264    }
1265    Int_t hh = TMath::Min(TMath::Nint(fraction * nPoints), Int_t(n));
1266    AliMathBase::EvaluateUni(nPoints, ddata, mean, lsigma, hh);
1267    if (sigma) *sigma = lsigma;
1268    delete [] ddata;
1269    return mean;
1270 }
1271
1272
1273 TH1F* AliTPCCalibViewer::SigmaCut(TH1F *const histogram, Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep, Bool_t pm) {
1274    //
1275    // 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
1276    // The data of the distribution Begin_Latex f(x, #mu, #sigma) End_Latex are given in 'histogram'
1277    // '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
1278    // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, Begin_Latex t #sigma End_Latex)
1279    // sigmaStep: the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
1280    // pm: Decide weather Begin_Latex t > 0 End_Latex (first case) or Begin_Latex t End_Latex arbitrary (secound case)
1281    // The actual work is done on the array.
1282    /* Begin_Latex 
1283          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    
1284          or      
1285          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 }
1286       End_Latex  
1287       begin_macro(source)
1288       {
1289          Float_t mean = 0;
1290          Float_t sigma = 1.5;
1291          Float_t sigmaMax = 4;
1292          gROOT->SetStyle("Plain");
1293          TH1F *distribution = new TH1F("Distribution1", "Distribution f(x, #mu, #sigma)", 1000,-5,5);
1294          TRandom rand(23);
1295          for (Int_t i = 0; i <50000;i++) distribution->Fill(rand.Gaus(mean, sigma));
1296          Float_t *ar = distribution->GetArray();
1297          
1298          TCanvas* macro_example_canvas = new TCanvas("macro_example_canvas_SigmaCut", "", 350, 350);
1299          macro_example_canvas->Divide(0,3);
1300          TVirtualPad *pad1 = macro_example_canvas->cd(1);
1301          pad1->SetGridy();
1302          pad1->SetGridx();
1303          distribution->Draw();
1304          TVirtualPad *pad2 = macro_example_canvas->cd(2);
1305          pad2->SetGridy();
1306          pad2->SetGridx();
1307          
1308          TH1F *shist = AliTPCCalibViewer::SigmaCut(distribution, mean, sigma, sigmaMax);
1309          shist->SetNameTitle("Cumulative","Cumulative S(t, #mu, #sigma)");
1310          shist->Draw();  
1311          TVirtualPad *pad3 = macro_example_canvas->cd(3);
1312          pad3->SetGridy();
1313          pad3->SetGridx();
1314          TH1F *shistPM = AliTPCCalibViewer::SigmaCut(distribution, mean, sigma, sigmaMax, -1, kTRUE);
1315          shistPM->Draw();   
1316          return macro_example_canvas;
1317       }  
1318       end_macro
1319    */ 
1320    
1321    Float_t *array = histogram->GetArray();
1322    Int_t    nbins = histogram->GetXaxis()->GetNbins();
1323    Float_t binLow = histogram->GetXaxis()->GetXmin();
1324    Float_t binUp  = histogram->GetXaxis()->GetXmax();
1325    return AliTPCCalibViewer::SigmaCut(nbins, array, mean, sigma, nbins, binLow, binUp, sigmaMax, sigmaStep, pm);
1326 }   
1327    
1328
1329 TH1F* AliTPCCalibViewer::SigmaCut(Int_t n, const 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){
1330    //
1331    // 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
1332    // The data of the distribution Begin_Latex f(x, #mu, #sigma) End_Latex are given in 'array', 'n' specifies the length of the array
1333    // '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
1334    // 'nbins': number of bins, 'binLow': first bin, 'binUp': last bin
1335    // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated (in units of sigma, Begin_Latex t #sigma End_Latex)
1336    // sigmaStep: the binsize of the generated histogram
1337    // Here the actual work is done.
1338    
1339    if (sigma == 0) return 0;
1340    Float_t binWidth = (binUp-binLow)/(nbins - 1);
1341    if (sigmaStep <= 0) sigmaStep = binWidth;
1342    Int_t kbins = (Int_t)(sigmaMax * sigma / sigmaStep) + 1; // + 1  due to overflow bin in histograms
1343    if (pm) kbins = 2 * (Int_t)(sigmaMax * sigma / sigmaStep) + 1;
1344    Float_t kbinLow = !pm ? 0 : -sigmaMax;
1345    Float_t kbinUp  = sigmaMax;
1346    TH1F *hist = new TH1F("sigmaCutHisto","Cumulative; Multiples of #sigma; Fraction of included data", kbins, kbinLow, kbinUp); 
1347    hist->SetDirectory(0);
1348    hist->Reset();
1349    
1350    // calculate normalization
1351    Double_t normalization = 0;
1352    for (Int_t i = 0; i <= n; i++) {
1353         normalization += array[i];
1354    }
1355    
1356    // given units: units from given histogram
1357    // sigma units: in units of sigma
1358    // iDelta: integrate in interval (mean +- iDelta), given units
1359    // x:      ofset from mean for integration, given units
1360    // hist:   needs 
1361    
1362 //    printf("nbins: %i, binLow: %f, binUp: %f \n", nbins, binLow, binUp);
1363    // fill histogram
1364    for (Float_t iDelta = 0; iDelta <= sigmaMax * sigma; iDelta += sigmaStep) {
1365       // integrate array
1366       Double_t valueP = array[GetBin(mean, nbins, binLow, binUp)];
1367       Double_t valueM = array[GetBin(mean-binWidth, nbins, binLow, binUp)];
1368       // add bin of mean value only once to the histogram
1369 //       printf("++ adding bins: ");
1370       for (Float_t x = binWidth; x <= iDelta; x += binWidth) {
1371          valueP += (mean + x <= binUp)  ? array[GetBin(mean + x, nbins, binLow, binUp)] : 0;
1372          valueM += (mean-binWidth - x >= binLow) ? array[GetBin(mean-binWidth - x, nbins, binLow, binUp)] : 0; 
1373 //          printf("%i, ", GetBin(mean + x, nbins, binLow, binUp));        
1374       }
1375 //       printf("\n");
1376       if (valueP / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail  +++ \n", valueP, normalization);
1377       if (valueP / normalization > 100) return hist;
1378       if (valueM / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail  +++ \n", valueM, normalization);
1379       if (valueM / normalization > 100) return hist;
1380       valueP = (valueP / normalization);
1381       valueM = (valueM / normalization);
1382       if (pm) {
1383          Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
1384          hist->SetBinContent(bin, valueP);
1385          bin = GetBin(-iDelta/sigma, kbins, kbinLow, kbinUp);
1386          hist->SetBinContent(bin, valueM);
1387       }
1388       else { // if (!pm)
1389          Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
1390          hist->SetBinContent(bin, valueP + valueM);
1391 //          printf("  first integration bin: %i, last integration bin in + direction: %i \n", GetBin(mean+binWidth, nbins, binLow, binUp), GetBin(iDelta, nbins, binLow, binUp));
1392 //          printf("  first integration bin: %i, last integration bin in - direction: %i \n", GetBin(mean+binWidth, nbins, binLow, binUp), GetBin(-iDelta, nbins, binLow, binUp));
1393 //          printf("  value: %f, normalization: %f, iDelta: %f, Bin: %i \n", valueP+valueM, normalization, iDelta, bin);
1394       }
1395    }
1396    //hist->SetMaximum(0.7);
1397    if (!pm) hist->SetMaximum(1.2);
1398    return hist;
1399 }
1400
1401
1402 TH1F* AliTPCCalibViewer::SigmaCut(Int_t n, const Double_t *array, Double_t mean, Double_t sigma, Int_t nbins, const Double_t *xbins, Double_t sigmaMax){
1403    // 
1404    // SigmaCut for variable binsize
1405    // NOT YET IMPLEMENTED !!!
1406    // 
1407    printf("SigmaCut with variable binsize, Not yet implemented\n");
1408    // avoid compiler warnings:
1409    n=n;
1410    mean=mean;
1411    sigma=sigma;
1412    nbins=nbins;
1413    sigmaMax=sigmaMax;
1414    array=array;
1415    xbins=xbins;
1416    
1417    return 0;
1418 }   
1419
1420
1421 TH1F* AliTPCCalibViewer::Integrate(TH1F *const histogram, Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep){
1422    //
1423    // 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"   
1424    // "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
1425    // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate 
1426    // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
1427    // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
1428    // The actual work is done on the array.
1429    /* Begin_Latex 
1430          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 }
1431       End_Latex  
1432       begin_macro(source)
1433       {
1434          Float_t mean = 0;
1435          Float_t sigma = 1.5;
1436          Float_t sigmaMax = 4;
1437          gROOT->SetStyle("Plain");
1438          TH1F *distribution = new TH1F("Distribution2", "Distribution f(x, #mu, #sigma)", 1000,-5,5);
1439          TRandom rand(23);
1440          for (Int_t i = 0; i <50000;i++) distribution->Fill(rand.Gaus(mean, sigma));
1441          Float_t *ar = distribution->GetArray();
1442          
1443          TCanvas* macro_example_canvas = new TCanvas("macro_example_canvas_Integrate", "", 350, 350);
1444          macro_example_canvas->Divide(0,2);
1445          TVirtualPad *pad1 = macro_example_canvas->cd(1);
1446          pad1->SetGridy();
1447          pad1->SetGridx();
1448          distribution->Draw();
1449          TVirtualPad *pad2 = macro_example_canvas->cd(2);
1450          pad2->SetGridy();
1451          pad2->SetGridx();
1452          TH1F *shist = AliTPCCalibViewer::Integrate(distribution, mean, sigma, sigmaMax);
1453          shist->SetNameTitle("Cumulative","Cumulative S(t, #mu, #sigma)");
1454          shist->Draw();  
1455          
1456          return macro_example_canvas_Integrate;
1457       }  
1458       end_macro
1459    */ 
1460
1461    
1462    Float_t *array = histogram->GetArray();
1463    Int_t    nbins = histogram->GetXaxis()->GetNbins();
1464    Float_t binLow = histogram->GetXaxis()->GetXmin();
1465    Float_t binUp  = histogram->GetXaxis()->GetXmax();
1466    return AliTPCCalibViewer::Integrate(nbins, array, nbins, binLow, binUp, mean, sigma, sigmaMax, sigmaStep);
1467 }   
1468
1469
1470 TH1F* AliTPCCalibViewer::Integrate(Int_t n, const Float_t *const array, Int_t nbins, Float_t binLow, Float_t binUp, Float_t mean, Float_t sigma, Float_t sigmaMax, Float_t sigmaStep){
1471    // 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"   
1472    // "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
1473    // sigmaMax: up to which sigma around the mean/median/LTM you want to integrate 
1474    // if "igma == 0" and "sigmaMax == 0" the whole histogram is integrated
1475    // "sigmaStep": the binsize of the generated histogram, -1 means, that the maximal reasonable stepsize is used
1476    // Here the actual work is done.
1477       
1478    Bool_t givenUnits = kTRUE;
1479    if (sigma != 0 && sigmaMax != 0) givenUnits = kFALSE;
1480    if (givenUnits) {
1481       sigma = 1;
1482       sigmaMax = (binUp - binLow) / 2.;
1483    }
1484    
1485    Float_t binWidth = (binUp-binLow)/(nbins - 1);
1486    if (sigmaStep <= 0) sigmaStep = binWidth;
1487    Int_t kbins =  (Int_t)(sigmaMax * sigma / sigmaStep) + 1;  // + 1  due to overflow bin in histograms
1488    Float_t kbinLow = givenUnits ? binLow : -sigmaMax;
1489    Float_t kbinUp  = givenUnits ? binUp  : sigmaMax;
1490    TH1F *hist = 0; 
1491    if (givenUnits)  hist = new TH1F("integratedHisto","Integrated Histogram; Given x; Fraction of included data", kbins, kbinLow, kbinUp); 
1492    if (!givenUnits) hist = new TH1F("integratedHisto","Integrated Histogram; Multiples of #sigma; Fraction of included data", kbins, kbinLow, kbinUp); 
1493    hist->SetDirectory(0);
1494    hist->Reset();
1495    
1496    // calculate normalization
1497  //  printf("calculating normalization, integrating from bin 1 to %i \n", n);
1498    Double_t normalization = 0;
1499    for (Int_t i = 1; i <= n; i++) {
1500         normalization += array[i];
1501    }
1502  //  printf("normalization: %f \n", normalization);
1503    
1504    // given units: units from given histogram
1505    // sigma units: in units of sigma
1506    // iDelta: integrate in interval (mean +- iDelta), given units
1507    // x:      ofset from mean for integration, given units
1508    // hist:   needs 
1509    
1510    // fill histogram
1511    for (Float_t iDelta = mean - sigmaMax * sigma; iDelta <= mean + sigmaMax * sigma; iDelta += sigmaStep) {
1512       // integrate array
1513       Double_t value = 0;
1514       for (Float_t x = mean - sigmaMax * sigma; x <= iDelta; x += binWidth) {
1515          value += (x <= binUp && x >= binLow)  ? array[GetBin(x, nbins, binLow, binUp)] : 0;
1516       }
1517       if (value / normalization > 100) printf("+++ Error, value to big: %f, normalization with %f will fail  +++ \n", value, normalization);
1518       if (value / normalization > 100) return hist;
1519       Int_t bin = GetBin(iDelta/sigma, kbins, kbinLow, kbinUp);
1520     //  printf("first integration bin: %i, last integration bin: %i \n", GetBin(mean - sigmaMax * sigma, nbins, binLow, binUp), GetBin(iDelta, nbins, binLow, binUp));
1521     //  printf("value: %f, normalization: %f, normalized value: %f, iDelta: %f, Bin: %i \n", value, normalization, value/normalization, iDelta, bin);
1522       value = (value / normalization);
1523       hist->SetBinContent(bin, value);
1524    }
1525    return hist;
1526 }
1527
1528
1529
1530
1531
1532 ////////////////////////
1533 // end of Array tools //
1534 ////////////////////////
1535
1536
1537
1538 //_____________________________________________________________________________
1539 AliTPCCalPad* AliTPCCalibViewer::GetCalPadOld(const char* desiredData, const char* cuts, const char* calPadName) const {
1540   //
1541   // creates a AliTPCCalPad out of the 'desiredData'
1542   // the functionality of EasyDraw1D is used
1543   // calPadName specifies the name of the created AliTPCCalPad
1544   //  - this takes a while -
1545   //
1546    TString drawStr(desiredData);
1547    drawStr.Append(":channel");
1548    drawStr.Append(fAbbreviation);
1549    AliTPCCalPad * createdCalPad = new AliTPCCalPad(calPadName, calPadName);
1550    Int_t entries = 0;
1551    for (Int_t sec = 0; sec < 72; sec++) {
1552      AliTPCCalROC * roc = createdCalPad->GetCalROC(sec);
1553       entries = EasyDraw1D(drawStr.Data(), (Int_t)sec, cuts, "goff");
1554       if (entries == -1) return 0;
1555       const Double_t *pchannel = fTree->GetV2();
1556       const Double_t *pvalue   = fTree->GetV1();
1557       for (Int_t i = 0; i < entries; i++) 
1558          roc->SetValue((UInt_t)(pchannel[i]), (Float_t)(pvalue[i]));
1559    }
1560    return createdCalPad;   
1561 }
1562
1563
1564 //_____________________________________________________________________________
1565 AliTPCCalPad* AliTPCCalibViewer::GetCalPad(const char* desiredData, const char* cuts, const char* calPadName) const {
1566   //
1567   // creates a AliTPCCalPad out of the 'desiredData'
1568   // the functionality of EasyDraw1D is used
1569   // calPadName specifies the name of the created AliTPCCalPad
1570   //  - this takes a while -
1571   //
1572    TString drawStr(desiredData);
1573    drawStr.Append(":channel.fElements:sector");
1574    AliTPCCalPad * createdCalPad = new AliTPCCalPad(calPadName, calPadName);
1575    //
1576    Int_t entries = fTree->Draw(drawStr, cuts,"goff");
1577    const Double_t *pvalue   = fTree->GetV1();
1578    const Double_t *pchannel = fTree->GetV2();
1579    const Double_t *psector  = fTree->GetV3();
1580
1581    for (Int_t ientry=0; ientry<entries; ientry++){
1582      Int_t sector= TMath::Nint(psector[ientry]);
1583      AliTPCCalROC * roc = createdCalPad->GetCalROC(sector);
1584      if (roc) roc->SetValue((UInt_t)(pchannel[ientry]), (Float_t)(pvalue[ientry]));
1585    }
1586
1587   //  for (Int_t sec = 0; sec < 72; sec++) {
1588 //      AliTPCCalROC * roc = createdCalPad->GetCalROC(sec);
1589 //       entries = EasyDraw1D(drawStr.Data(), (Int_t)sec, cuts, "goff");
1590 //       if (entries == -1) return 0;
1591 //       for (Int_t i = 0; i < entries; i++) 
1592 //          roc->SetValue((UInt_t)(pchannel[i]), (Float_t)(pvalue[i]));
1593 //    }
1594    return createdCalPad;   
1595 }
1596
1597 //_____________________________________________________________________________
1598 AliTPCCalROC* AliTPCCalibViewer::GetCalROC(const char* desiredData, UInt_t sector, const char* cuts) const {
1599   //
1600   // creates a AliTPCCalROC out of the desiredData
1601   // the functionality of EasyDraw1D is used
1602   // sector specifies the sector of the created AliTPCCalROC
1603   //
1604    TString drawStr(desiredData);
1605    drawStr.Append(":channel");
1606    drawStr.Append(fAbbreviation);
1607    Int_t entries = EasyDraw1D(drawStr.Data(), (Int_t)sector, cuts, "goff");
1608    if (entries == -1) return 0;
1609    AliTPCCalROC * createdROC = new AliTPCCalROC(sector);
1610    for (Int_t i = 0; i < entries; i++) 
1611       createdROC->SetValue((UInt_t)(fTree->GetV2()[i]), fTree->GetV1()[i]);
1612    return createdROC;
1613 }
1614
1615
1616 TObjArray* AliTPCCalibViewer::GetListOfVariables(Bool_t printList) {
1617   //
1618   // scan the tree  - produces a list of available variables in the tree
1619   // printList: print the list to the screen, after the scan is done
1620   //
1621    TObjArray* arr = new TObjArray();
1622    TObjString* str = 0;
1623    if (!fTree) return 0;
1624    Int_t nentries = fTree->GetListOfBranches()->GetEntries();
1625    for (Int_t i = 0; i < nentries; i++) {
1626       str = new TObjString(fTree->GetListOfBranches()->At(i)->GetName());
1627       str->String().ReplaceAll("_Median", "");
1628       str->String().ReplaceAll("_Mean", "");
1629       str->String().ReplaceAll("_RMS", "");
1630       str->String().ReplaceAll("_LTM", "");
1631       str->String().ReplaceAll("_OutlierCutted", "");
1632       str->String().ReplaceAll(".", "");
1633       if (!arr->FindObject(str) && 
1634           !(str->String() == "channel" || str->String() == "gx" || str->String() == "gy" || 
1635             str->String() == "lx" || str->String() == "ly" || str->String() == "pad" || 
1636             str->String() == "row" || str->String() == "rpad" || str->String() == "sector"  ))
1637          arr->Add(str);
1638    }
1639    
1640    // loop over all friends (if there are some) and add them to the list
1641    if (fTree->GetListOfFriends()) {
1642       for (Int_t ifriend = 0; ifriend < fTree->GetListOfFriends()->GetEntries(); ifriend++){
1643          // printf("iterating through friendlist, currently at %i\n", ifriend);
1644          // printf("working with %s\n", fTree->GetListOfFriends()->At(ifriend)->ClassName());
1645          if (TString(fTree->GetListOfFriends()->At(ifriend)->ClassName()) != "TFriendElement") continue; // no friendElement found
1646          TFriendElement *friendElement = (TFriendElement*)fTree->GetListOfFriends()->At(ifriend);
1647          if (friendElement->GetTree() == 0) continue; // no tree found in friendElement
1648          // printf("friend found \n");
1649          for (Int_t i = 0; i < friendElement->GetTree()->GetListOfBranches()->GetEntries(); i++) {
1650             // printf("iterating through friendelement entries, currently at %i\n", i);
1651             str = new TObjString(friendElement->GetTree()->GetListOfBranches()->At(i)->GetName());
1652             str->String().ReplaceAll("_Median", "");
1653             str->String().ReplaceAll("_Mean", "");
1654             str->String().ReplaceAll("_RMS", "");
1655             str->String().ReplaceAll("_LTM", "");
1656             str->String().ReplaceAll("_OutlierCutted", "");
1657             str->String().ReplaceAll(".", "");
1658             if (!(str->String() == "channel" || str->String() == "gx" || str->String() == "gy" || 
1659                   str->String() == "lx" || str->String() == "ly" || str->String() == "pad" || 
1660                   str->String() == "row" || str->String() == "rpad" || str->String() == "sector"  )){
1661                // insert "<friendName>." at the beginning: (<friendName> is per default "R")
1662                str->String().Insert(0, ".");
1663                str->String().Insert(0, friendElement->GetName());
1664                if (!arr->FindObject(str)) arr->Add(str);
1665                // printf("added string %s \n", str->String().Data());
1666             }
1667          }
1668       }
1669    } // if (fTree->GetListOfFriends())
1670    
1671    arr->Sort();
1672 //   ((TFriendElement*)gui->GetViewer()->GetTree()->GetListOfFriends()->At(0))->GetTree()->GetListOfBranches()->At(0)->GetName()
1673 // ((TFriendElement*)gui->GetViewer()->GetTree()->GetListOfFriends()->At(0))->GetTree()->GetListOfBranches()
1674
1675
1676    if (printList) {
1677       TIterator* iter = arr->MakeIterator();
1678       iter->Reset();
1679       TObjString* currentStr = 0;
1680       while ( (currentStr = (TObjString*)(iter->Next())) ) {
1681          std::cout << currentStr->GetString().Data() << std::endl;
1682       }
1683       delete iter;
1684    }
1685    return arr;
1686 }
1687
1688
1689 TObjArray* AliTPCCalibViewer::GetListOfNormalizationVariables(Bool_t printList) const{
1690   //
1691   // produces a list of available variables for normalization in the tree
1692   // printList: print the list to the screen, after the scan is done
1693   //
1694    TObjArray* arr = new TObjArray();
1695    arr->Add(new TObjString("_Mean"));
1696    arr->Add(new TObjString("_Mean_OutlierCutted"));
1697    arr->Add(new TObjString("_Median"));
1698    arr->Add(new TObjString("_Median_OutlierCutted"));
1699    arr->Add(new TObjString("_LTM"));
1700    arr->Add(new TObjString("_LTM_OutlierCutted"));
1701    arr->Add(new TObjString(Form("LFitIntern_4_8%s", fAppendString.Data())));
1702    arr->Add(new TObjString(Form("GFitIntern_Lin%s", fAppendString.Data())));
1703    arr->Add(new TObjString(Form("GFitIntern_Par%s", fAppendString.Data())));
1704    arr->Add(new TObjString("FitLinLocal"));
1705    arr->Add(new TObjString("FitLinGlobal"));
1706    arr->Add(new TObjString("FitParLocal"));
1707    arr->Add(new TObjString("FitParGlobal"));
1708
1709    if (printList) {
1710       TIterator* iter = arr->MakeIterator();
1711       iter->Reset();
1712       TObjString* currentStr = 0;
1713       while ((currentStr = (TObjString*)(iter->Next()))) {
1714          std::cout << currentStr->GetString().Data() << std::endl;
1715       }
1716       delete iter;
1717    }
1718    return arr;
1719 }
1720
1721
1722 TFriendElement* AliTPCCalibViewer::AddReferenceTree(const char* filename, const char* treename, const char* refname){
1723   //
1724   // add a reference tree to the current tree
1725   // by default the treename is 'calPads' and the reference treename is 'R'
1726   //
1727    TFile *file = new TFile(filename);
1728    fListOfObjectsToBeDeleted->Add(file);
1729    TTree * tree = (TTree*)file->Get(treename);
1730    return AddFriend(tree, refname);
1731 }
1732
1733
1734 TObjArray* AliTPCCalibViewer::GetArrayOfCalPads(){
1735   //
1736   // Returns a TObjArray with all AliTPCCalPads that are stored in the tree
1737   //  - this takes a while - 
1738   //
1739    TObjArray *listOfCalPads = GetListOfVariables();
1740    TObjArray *calPadsArray = new TObjArray();
1741    Int_t numberOfCalPads = listOfCalPads->GetEntries();
1742    for (Int_t i = 0; i < numberOfCalPads; i++) {
1743       std::cout << "Creating calPad " << (i+1) << " of " << numberOfCalPads << "\r" << std::flush;
1744       char* calPadName = (char*)((TObjString*)(listOfCalPads->At(i)))->GetString().Data();
1745       TString drawCommand = ((TObjString*)(listOfCalPads->At(i)))->GetString();
1746       drawCommand.Append(fAbbreviation.Data());
1747       AliTPCCalPad* calPad = GetCalPad(drawCommand.Data(), "", calPadName); 
1748       calPadsArray->Add(calPad); 
1749    }
1750    std::cout << std::endl;
1751    listOfCalPads->Delete();
1752    delete listOfCalPads;
1753    return calPadsArray;
1754 }
1755
1756
1757 TString* AliTPCCalibViewer::Fit(const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, TVectorD &fitParam, TMatrixD &covMatrix){
1758    //
1759    // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts
1760    // returns chi2, fitParam and covMatrix
1761    // returns TString with fitted formula
1762    //
1763    
1764    TString formulaStr(formula); 
1765    TString drawStr(drawCommand);
1766    TString cutStr(cuts);
1767    
1768    // abbreviations:
1769    drawStr.ReplaceAll(fAbbreviation, fAppendString);
1770    cutStr.ReplaceAll(fAbbreviation, fAppendString);
1771    formulaStr.ReplaceAll(fAbbreviation, fAppendString);
1772    
1773    formulaStr.ReplaceAll("++", fAbbreviation);
1774    TObjArray* formulaTokens = formulaStr.Tokenize(fAbbreviation.Data()); 
1775    Int_t dim = formulaTokens->GetEntriesFast();
1776    
1777    fitParam.ResizeTo(dim);
1778    covMatrix.ResizeTo(dim,dim);
1779    
1780    TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim));
1781    fitter->StoreData(kTRUE);   
1782    fitter->ClearPoints();
1783    
1784    Int_t entries = Draw(drawStr.Data(), cutStr.Data(), "goff");
1785    if (entries == -1) return new TString("An ERROR has occured during fitting!");
1786    Double_t **values = new Double_t*[dim+1] ; 
1787    
1788    for (Int_t i = 0; i < dim + 1; i++){
1789       Int_t centries = 0;
1790       if (i < dim) centries = fTree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff");
1791       else  centries = fTree->Draw(drawStr.Data(), cutStr.Data(), "goff");
1792       
1793       if (entries != centries) return new TString("An ERROR has occured during fitting!");
1794       values[i] = new Double_t[entries];
1795       memcpy(values[i],  fTree->GetV1(), entries*sizeof(Double_t)); 
1796    }
1797    
1798    // add points to the fitter
1799    for (Int_t i = 0; i < entries; i++){
1800       Double_t x[1000];
1801       for (Int_t j=0; j<dim;j++) x[j]=values[j][i];
1802       fitter->AddPoint(x, values[dim][i], 1);
1803    }
1804
1805    fitter->Eval();
1806    fitter->GetParameters(fitParam);
1807    fitter->GetCovarianceMatrix(covMatrix);
1808    chi2 = fitter->GetChisquare();
1809    chi2 = chi2;
1810    
1811    TString *preturnFormula = new TString(Form("( %e+",fitParam[0])), &returnFormula = *preturnFormula;
1812    
1813    for (Int_t iparam = 0; iparam < dim; iparam++) {
1814      returnFormula.Append(Form("%s*(%e)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1]));
1815      if (iparam < dim-1) returnFormula.Append("+");
1816    }
1817    returnFormula.Append(" )");
1818    delete formulaTokens;
1819    delete fitter;
1820    delete[] values;
1821    return preturnFormula;
1822 }
1823
1824
1825 void AliTPCCalibViewer::MakeTreeWithObjects(const char *fileName, const TObjArray *const array, const char * mapFileName) {
1826   //
1827   // Write tree with all available information
1828   // im mapFileName is speciefied, the Map information are also written to the tree
1829   // AliTPCCalPad-Objects are written directly to the tree, so that they can be accessd later on
1830   // (does not work!!!)
1831   //
1832    AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
1833
1834    TObjArray* mapIROCs = 0;
1835    TObjArray* mapOROCs = 0;
1836    TVectorF *mapIROCArray = 0;
1837    TVectorF *mapOROCArray = 0;
1838    Int_t mapEntries = 0;
1839    TString* mapNames = 0;
1840    
1841    if (mapFileName) {
1842       TFile mapFile(mapFileName, "read");
1843       
1844       TList* listOfROCs = mapFile.GetListOfKeys();
1845       mapEntries = listOfROCs->GetEntries()/2;
1846       mapIROCs = new TObjArray(mapEntries*2);
1847       mapOROCs = new TObjArray(mapEntries*2);
1848       mapIROCArray = new TVectorF[mapEntries];
1849       mapOROCArray = new TVectorF[mapEntries];
1850       
1851       mapNames = new TString[mapEntries];
1852       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1853          TString rocName(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
1854          rocName.Remove(rocName.Length()-4, 4);
1855          mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((rocName + "IROC").Data()), ivalue);
1856          mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((rocName + "OROC").Data()), ivalue);
1857          mapNames[ivalue].Append(rocName);
1858       }
1859       
1860       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1861          mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
1862          mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
1863       
1864          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
1865             (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
1866          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
1867             (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
1868       }
1869
1870    } //  if (mapFileName)
1871   
1872    TTreeSRedirector cstream(fileName);
1873    Int_t arrayEntries = array->GetEntries();
1874    
1875    // Read names of AliTPCCalPads and save them in names[]
1876    TString* names = new TString[arrayEntries];
1877    for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
1878       names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
1879
1880    for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
1881       
1882       TVectorF *vectorArray = new TVectorF[arrayEntries];
1883       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
1884          vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
1885             
1886       
1887       //
1888       // fill vectors of variable per pad
1889       //
1890       TVectorF *posArray = new TVectorF[8];
1891       for (Int_t ivalue = 0; ivalue < 8; ivalue++)
1892          posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
1893
1894       Float_t posG[3] = {0};
1895       Float_t posL[3] = {0};
1896       Int_t ichannel = 0;
1897       for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
1898          for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
1899             tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
1900             tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
1901             posArray[0][ichannel] = irow;
1902             posArray[1][ichannel] = ipad;
1903             posArray[2][ichannel] = posL[0];
1904             posArray[3][ichannel] = posL[1];
1905             posArray[4][ichannel] = posG[0];
1906             posArray[5][ichannel] = posG[1];
1907             posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
1908             posArray[7][ichannel] = ichannel;
1909             
1910             // loop over array containing AliTPCCalPads
1911             for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
1912                AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
1913                AliTPCCalROC* calROC = calPad->GetCalROC(isector);
1914                if (calROC)
1915                   (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
1916                else
1917                   (vectorArray[ivalue])[ichannel] = 0;
1918             }
1919             ichannel++;
1920          }
1921       }
1922       AliTPCCalROC dummyROC(0);
1923       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
1924          AliTPCCalROC *roc = ((AliTPCCalPad*)array->At(ivalue))->GetCalROC(isector);
1925          if (!roc) roc = &dummyROC;
1926          cstream << "calPads" <<
1927             (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
1928          cstream << "calPads" << 
1929             (Char_t*)((names[ivalue] + "Pad.=").Data()) << roc;
1930       }
1931
1932       if (mapFileName) {
1933          for  (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
1934             if (isector < 36)
1935                cstream << "calPads" <<
1936                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
1937             else
1938                cstream << "calPads" <<
1939                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
1940          }
1941       }
1942       
1943       cstream << "calPads" <<
1944          "sector=" << isector;
1945
1946       cstream << "calPads" <<
1947          "row.=" << &posArray[0] <<
1948          "pad.=" << &posArray[1] <<
1949          "lx.=" << &posArray[2] <<
1950          "ly.=" << &posArray[3] <<
1951          "gx.=" << &posArray[4] <<
1952          "gy.=" << &posArray[5] <<
1953          "rpad.=" << &posArray[6] <<
1954          "channel.=" << &posArray[7];
1955
1956       cstream << "calPads" <<
1957          "\n";
1958
1959       delete[] posArray;
1960       delete[] vectorArray;
1961    } //for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++)
1962
1963    delete[] names;
1964    if (mapFileName) {
1965       delete mapIROCs;
1966       delete mapOROCs;
1967       delete[] mapIROCArray;
1968       delete[] mapOROCArray;
1969       delete[] mapNames;
1970    }
1971 }
1972
1973
1974 void AliTPCCalibViewer::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad *const outlierPad, Float_t ltmFraction) {
1975   //
1976   // Write a tree with all available information
1977   // if mapFileName is speciefied, the Map information are also written to the tree
1978   // pads specified in outlierPad are not used for calculating statistics
1979   // The following statistical information on the basis of a ROC are calculated: 
1980   // "_Median", "_Mean", "_LTM", "_RMS_LTM"
1981   // "_Median_OutlierCutted", "_Mean_OutlierCutted", "_RMS_OutlierCutted", "_LTM_OutlierCutted", "_RMS_LTM_OutlierCutted"
1982   // The following position variables are available:
1983   // "row", "pad", "lx", "ly", "gx", "gy", "rpad", "channel"
1984   // 
1985   // The tree out of this function is the basis for the AliTPCCalibViewer and the AliTPCCalibViewerGUI.
1986    
1987   AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
1988
1989   TObjArray* mapIROCs = 0;
1990   TObjArray* mapOROCs = 0;
1991   TVectorF *mapIROCArray = 0;
1992   TVectorF *mapOROCArray = 0;
1993   Int_t mapEntries = 0;
1994   TString* mapNames = 0;
1995   
1996   if (mapFileName) {
1997     TFile mapFile(mapFileName, "read");
1998     
1999     TList* listOfROCs = mapFile.GetListOfKeys();
2000     mapEntries = listOfROCs->GetEntries()/2;
2001     mapIROCs = new TObjArray(mapEntries*2);
2002     mapOROCs = new TObjArray(mapEntries*2);
2003     mapIROCArray = new TVectorF[mapEntries];
2004     mapOROCArray = new TVectorF[mapEntries];
2005     
2006     mapNames = new TString[mapEntries];
2007     for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
2008       TString rocName(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
2009       rocName.Remove(rocName.Length()-4, 4);
2010       mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((rocName + "IROC").Data()), ivalue);
2011       mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((rocName + "OROC").Data()), ivalue);
2012       mapNames[ivalue].Append(rocName);
2013     }
2014     
2015     for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
2016       mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
2017       mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
2018       
2019       for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
2020         (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
2021       for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
2022         (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
2023     }
2024     
2025   } //  if (mapFileName)
2026   
2027   TTreeSRedirector cstream(fileName);
2028   Int_t arrayEntries = 0;
2029   if (array) arrayEntries = array->GetEntries();
2030   
2031   TString* names = new TString[arrayEntries];
2032   for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
2033     names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
2034   
2035   for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
2036       //
2037       // get statistic for given sector
2038       //
2039     TVectorF median(arrayEntries);
2040     TVectorF mean(arrayEntries);
2041     TVectorF rms(arrayEntries);
2042     TVectorF ltm(arrayEntries);
2043     TVectorF ltmrms(arrayEntries);
2044     TVectorF medianWithOut(arrayEntries);
2045     TVectorF meanWithOut(arrayEntries);
2046     TVectorF rmsWithOut(arrayEntries);
2047     TVectorF ltmWithOut(arrayEntries);
2048     TVectorF ltmrmsWithOut(arrayEntries);
2049     
2050     TVectorF *vectorArray = new TVectorF[arrayEntries];
2051     for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++){
2052       vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
2053       vectorArray[ivalue].SetUniqueID(array->UncheckedAt(ivalue)->GetUniqueID());
2054 //       printf("UniqueID: %d\n",vectorArray[ivalue].GetUniqueID());
2055     }
2056     
2057     for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
2058       AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
2059       AliTPCCalROC* calROC = calPad->GetCalROC(isector);
2060       AliTPCCalROC* outlierROC = 0;
2061       if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
2062       if (calROC) {
2063         median[ivalue] = calROC->GetMedian();
2064         mean[ivalue] = calROC->GetMean();
2065         rms[ivalue] = calROC->GetRMS();
2066         Double_t ltmrmsValue = 0;
2067         ltm[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction);
2068         ltmrms[ivalue] = ltmrmsValue;
2069         if (outlierROC) {
2070           medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
2071           meanWithOut[ivalue] = calROC->GetMean(outlierROC);
2072           rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
2073           ltmrmsValue = 0;
2074           ltmWithOut[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction, outlierROC);
2075           ltmrmsWithOut[ivalue] = ltmrmsValue;
2076         }
2077       }
2078       else {
2079         median[ivalue] = 0.;
2080         mean[ivalue] = 0.;
2081         rms[ivalue] = 0.;
2082         ltm[ivalue] = 0.;
2083         ltmrms[ivalue] = 0.;
2084         medianWithOut[ivalue] = 0.;
2085         meanWithOut[ivalue] = 0.;
2086         rmsWithOut[ivalue] = 0.;
2087         ltmWithOut[ivalue] = 0.;
2088         ltmrmsWithOut[ivalue] = 0.;
2089       }
2090     }
2091     
2092       //
2093       // fill vectors of variable per pad
2094       //
2095     TVectorF *posArray = new TVectorF[8];
2096     for (Int_t ivalue = 0; ivalue < 8; ivalue++)
2097       posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
2098     
2099     Float_t posG[3] = {0};
2100     Float_t posL[3] = {0};
2101     Int_t ichannel = 0;
2102     for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
2103       for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
2104         tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
2105         tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
2106         posArray[0][ichannel] = irow;
2107         posArray[1][ichannel] = ipad;
2108         posArray[2][ichannel] = posL[0];
2109         posArray[3][ichannel] = posL[1];
2110         posArray[4][ichannel] = posG[0];
2111         posArray[5][ichannel] = posG[1];
2112         posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
2113         posArray[7][ichannel] = ichannel;
2114         
2115             // loop over array containing AliTPCCalPads
2116         for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
2117           AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
2118           AliTPCCalROC* calROC = calPad->GetCalROC(isector);
2119                 if (calROC)
2120                   (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
2121           else
2122             (vectorArray[ivalue])[ichannel] = 0;
2123         }
2124         ichannel++;
2125       }
2126     }
2127     
2128     cstream << "calPads" <<
2129       "sector=" << isector;
2130     
2131     for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
2132       cstream << "calPads" <<
2133         (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
2134         (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
2135         (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
2136         (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
2137         (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
2138       if (outlierPad) {
2139         cstream << "calPads" <<
2140           (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
2141           (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
2142           (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
2143           (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
2144           (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
2145       }
2146         //timestamp and run, if given in title
2147 /*      TString title(((AliTPCCalPad*) array->At(ivalue))->GetTitle());
2148       TObjArray *arrtitle=title.Tokenize(",");
2149       Int_t run=-1;
2150       UInt_t time=0;
2151       TIter next(arrtitle);
2152       TObject *o=0;
2153       while ( (o=next()) ){
2154         TString &entry=((TObjString*)o)->GetString();
2155         entry.Remove(TString::kBoth,' ');
2156         entry.Remove(TString::kBoth,'\t');
2157         if (entry.BeginsWith("Run:")) {
2158           run=entry(4,entry.Length()).Atoi();
2159         } else if (entry.BeginsWith("Time:")) {
2160           time=entry(6,entry.Length()).Atoi();
2161         }
2162       }
2163       delete arrtitle;*/
2164       
2165     }
2166     
2167     for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
2168       cstream << "calPads" <<
2169         (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
2170     }
2171     
2172     if (mapFileName) {
2173           for  (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
2174             if (isector < 36)
2175               cstream << "calPads" <<
2176               (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
2177             else
2178                 cstream << "calPads" <<
2179               (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
2180          }
2181     }
2182
2183     cstream << "calPads" <<
2184       "row.=" << &posArray[0] <<
2185       "pad.=" << &posArray[1] <<
2186       "lx.=" << &posArray[2] <<
2187       "ly.=" << &posArray[3] <<
2188       "gx.=" << &posArray[4] <<
2189       "gy.=" << &posArray[5] <<
2190       "rpad.=" << &posArray[6] <<
2191       "channel.=" << &posArray[7];
2192     
2193     cstream << "calPads" <<
2194       "\n";
2195     
2196     delete[] posArray;
2197     delete[] vectorArray;
2198   }
2199   
2200   
2201   delete[] names;
2202   if (mapFileName) {
2203     delete mapIROCs;
2204     delete mapOROCs;
2205     delete[] mapIROCArray;
2206     delete[] mapOROCArray;
2207     delete[] mapNames;
2208   }
2209 }
2210
2211
2212 void AliTPCCalibViewer::MakeTree(const char *outPutFileName, const Char_t *inputFileName, AliTPCCalPad *outlierPad, Float_t ltmFraction, const char *mapFileName ){
2213    // 
2214    // Function to create a calibration Tree with all available information.
2215    // See also documentation to MakeTree   
2216    // the file "inputFileName" must be in the following format (see also CreateObjectList):
2217    // (each colum separated by tabs, "dependingOnType" can have 2 or 3 colums)
2218    // 
2219    // type      path    dependingOnType
2220    // 
2221    // type == "CE":
2222    // dependingOnType = CETmean CEQmean CETrms
2223    // 
2224    // type == "Pulser":
2225    // dependingOnType = PulserTmean     PulsterQmean    PulserTrms
2226    // 
2227    // type == "Pedestals":
2228    // dependingOnType = Pedestals       Noise
2229    // 
2230    // type == "CalPad":
2231    // dependingOnType = NameInFile      NameToWriteToFile
2232    // 
2233    // 
2234    TObjArray objArray;
2235    CreateObjectList(inputFileName, &objArray);
2236    MakeTree(outPutFileName, &objArray, mapFileName, outlierPad, ltmFraction);   
2237 }
2238
2239
2240 void AliTPCCalibViewer::CreateObjectList(const Char_t *filename, TObjArray *calibObjects){
2241    // 
2242    // Function to create a TObjArray out of a given file
2243    // the file must be in the following format:
2244    // (each colum separated by tabs, "dependingOnType" can have 2 or 3 colums)
2245    // 
2246    // 
2247    // type      path    dependingOnType
2248    // 
2249    // type == "CE":
2250    // dependingOnType = CETmean CEQmean CETrms
2251    // 
2252    // type == "Pulser":
2253    // dependingOnType = PulserTmean     PulsterQmean    PulserTrms
2254    // 
2255    // type == "Pedestals":
2256    // dependingOnType = Pedestals       Noise
2257    // 
2258    // type == "CalPad":
2259    // dependingOnType = NameInFile      NameToWriteToFile
2260    // 
2261    // 
2262    // 
2263    if ( calibObjects == 0x0 ) return;
2264    ifstream in;
2265    in.open(filename);
2266    if ( !in.is_open() ){
2267       fprintf(stderr,"Error: cannot open list file '%s'", filename);
2268       return;
2269    }
2270    
2271    AliTPCCalPad *calPad=0x0;
2272    
2273    TString sFile;
2274    sFile.ReadFile(in);
2275    in.close();
2276    
2277    TObjArray *arrFileLine = sFile.Tokenize("\n");
2278    TIter nextLine(arrFileLine);
2279    
2280    TObjString *sObjLine = 0x0;
2281    while ( (sObjLine = (TObjString*)nextLine()) ){
2282       TString sLine(sObjLine->GetString());
2283       
2284       TObjArray *arrCol = sLine.Tokenize("\t");
2285       Int_t nCols = arrCol->GetEntriesFast();
2286       
2287       TObjString *sObjType     = (TObjString*)(arrCol->At(0));
2288       TObjString *sObjFileName = (TObjString*)(arrCol->At(1));
2289       TObjString *sObjName = 0x0;
2290       
2291       if ( !sObjType || !sObjFileName ) continue;
2292       TString sType(sObjType->GetString());
2293       TString sFileName(sObjFileName->GetString());
2294       printf("Type %s, opening %s \n", sType.Data(), sFileName.Data());
2295       TFile *fIn = TFile::Open(sFileName);
2296       if ( !fIn ){
2297          fprintf(stderr,"File not found: '%s'", sFileName.Data());
2298          continue;
2299       }
2300       
2301       if ( sType == "CE" ){  // next three colums are the names for CETmean, CEQmean and CETrms
2302          AliTPCCalibCE *ce = (AliTPCCalibCE*)fIn->Get("AliTPCCalibCE");
2303          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadT0());         
2304          if (nCols > 2) {  // check, if the name is provided
2305             sObjName = (TObjString*)(arrCol->At(2));
2306             calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2307          }
2308          else calPad->SetNameTitle("CETmean","CETmean");
2309          calibObjects->Add(calPad);
2310          
2311          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadQ());         
2312          if (nCols > 3) {  // check, if the name is provided
2313             sObjName = (TObjString*)(arrCol->At(3));
2314             calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2315          }
2316          else calPad->SetNameTitle("CEQmean","CEQmean");
2317          calibObjects->Add(calPad);        
2318          
2319          calPad = new AliTPCCalPad((TObjArray*)ce->GetCalPadRMS());
2320          if (nCols > 4) {  // check, if the name is provided
2321             sObjName = (TObjString*)(arrCol->At(4));
2322             calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2323          }
2324          else calPad->SetNameTitle("CETrms","CETrms");
2325          calibObjects->Add(calPad);         
2326                   
2327       } else if ( sType == "Pulser") {
2328          AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fIn->Get("AliTPCCalibPulser");
2329          
2330          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadT0());         
2331          if (nCols > 2) {
2332             sObjName = (TObjString*)(arrCol->At(2));
2333             calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2334          }
2335          else calPad->SetNameTitle("PulserTmean","PulserTmean");
2336          calibObjects->Add(calPad);
2337          
2338          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadQ());         
2339          if (nCols > 3) {
2340             sObjName = (TObjString*)(arrCol->At(3));
2341             calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2342          }
2343          else calPad->SetNameTitle("PulserQmean","PulserQmean");
2344          calibObjects->Add(calPad);        
2345          
2346          calPad = new AliTPCCalPad((TObjArray*)sig->GetCalPadRMS());
2347          if (nCols > 4) {
2348             sObjName = (TObjString*)(arrCol->At(4));
2349             calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2350          }
2351          else calPad->SetNameTitle("PulserTrms","PulserTrms");
2352          calibObjects->Add(calPad);         
2353       
2354       } else if ( sType == "Pedestals") {
2355          AliTPCCalibPedestal *ped = (AliTPCCalibPedestal*)fIn->Get("AliTPCCalibPedestal");
2356          
2357          calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadPedestal());         
2358          if (nCols > 2) {
2359             sObjName = (TObjString*)(arrCol->At(2));
2360             calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2361          }
2362          else calPad->SetNameTitle("Pedestals","Pedestals");
2363          calibObjects->Add(calPad);
2364          
2365          calPad = new AliTPCCalPad((TObjArray*)ped->GetCalPadRMS());         
2366          if (nCols > 3) {
2367             sObjName = (TObjString*)(arrCol->At(3));
2368             calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2369          }
2370          else calPad->SetNameTitle("Noise","Noise");
2371          calibObjects->Add(calPad);        
2372      
2373       } else if ( sType == "CalPad") {
2374          if (nCols > 2) sObjName = (TObjString*)(arrCol->At(2));
2375          else continue;
2376          calPad = new AliTPCCalPad(*(AliTPCCalPad*)fIn->Get(sObjName->GetString().Data()));
2377          if (nCols > 3) {
2378             sObjName = (TObjString*)(arrCol->At(3));
2379             calPad->SetNameTitle(sObjName->GetString().Data(), sObjName->GetString().Data());
2380          }
2381          calibObjects->Add(calPad);
2382       } else {
2383          fprintf(stderr,"Undefined Type: '%s'",sType.Data());
2384       }
2385       delete fIn;
2386    }
2387 }
2388
2389
2390