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