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