]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCalibViewer.cxx
f71a4db184b6d28b9afa138bf58e405375589c6b
[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
23 //
24 // ROOT includes 
25 //
26 #include <iostream>
27 #include <TString.h>
28 #include <TRandom.h>
29 #include <TLegend.h>
30 #include <TLine.h>
31 #include <TCanvas.h>
32 #include <TROOT.h>
33 #include <TStyle.h>
34 #include <TH1F.h>
35 #include <THashTable.h>
36 #include <TObjString.h>
37 #include "TTreeStream.h"
38 #include "TFile.h"
39 #include "TKey.h"
40
41
42 //
43 // AliRoot includes
44 //
45 #include "AliTPCCalibViewer.h"
46
47 ClassImp(AliTPCCalibViewer)
48
49 AliTPCCalibViewer::AliTPCCalibViewer()
50                   :TObject(),
51                    fTree(0),
52                    fFile(0),
53                    fListOfObjectsToBeDeleted(0)
54 {
55   //
56   // Default constructor
57   //
58
59 }
60
61 //_____________________________________________________________________________
62 AliTPCCalibViewer::AliTPCCalibViewer(const AliTPCCalibViewer &c)
63                   :TObject(c),
64                    fTree(0),
65                    fFile(0),
66                    fListOfObjectsToBeDeleted(0)
67 {
68   //
69   // dummy AliTPCCalibViewer copy constructor
70   // not yet working!!!
71   //
72   fTree = c.fTree;
73   //fFile = new TFile(*(c.fFile));
74   fListOfObjectsToBeDeleted = c.fListOfObjectsToBeDeleted;
75 }
76
77 //_____________________________________________________________________________
78 AliTPCCalibViewer::AliTPCCalibViewer(TTree* tree)
79                   :TObject(),
80                    fTree(0),
81                    fFile(0),
82                    fListOfObjectsToBeDeleted(0)
83 {
84   //
85   // Constructor that initializes the calibration viewer
86   //
87   fTree = tree;
88   fListOfObjectsToBeDeleted = new TObjArray();
89 }
90
91 //_____________________________________________________________________________
92 AliTPCCalibViewer::AliTPCCalibViewer(char* fileName, char* treeName)
93                   :TObject(),
94                    fTree(0),
95                    fFile(0),
96                    fListOfObjectsToBeDeleted(0)
97 {
98    //
99    // Constructor to initialize the calibration viewer
100    // the file 'fileName' contains the tree 'treeName'
101    //
102    fFile = new TFile(fileName, "read");
103    fTree = (TTree*) fFile->Get(treeName);
104    fListOfObjectsToBeDeleted = new TObjArray();
105 }
106                    
107 //____________________________________________________________________________
108 AliTPCCalibViewer & AliTPCCalibViewer::operator =(const AliTPCCalibViewer & param)
109 {
110    //
111    // assignment operator - dummy
112    // not yet working!!!
113    //
114    fTree = param.fTree;
115    //fFile = new TFile(*(param.fFile));
116    fListOfObjectsToBeDeleted = param.fListOfObjectsToBeDeleted;
117    return (*this);
118 }
119
120 //_____________________________________________________________________________
121 AliTPCCalibViewer::~AliTPCCalibViewer()
122 {
123    //
124    // AliTPCCalibViewer destructor
125    // all objects will be deleted, the file will be closed, the pictures will disapear
126    //
127    /*if (fTree) {
128       delete fTree;
129       fTree = 0;
130    }*/
131    if (fFile) {
132       fFile->Close();
133       fFile = 0;
134    }
135
136    for (Int_t i = fListOfObjectsToBeDeleted->GetEntriesFast()-1; i >= 0; i--) {
137       //cout << "Index " << i << " trying to delete the following object: " << fListOfObjectsToBeDeleted->At(i)->GetName() << "..."<< endl;
138       delete fListOfObjectsToBeDeleted->At(i);
139    }
140    delete fListOfObjectsToBeDeleted;
141 }
142
143 //_____________________________________________________________________________
144 Int_t AliTPCCalibViewer::EasyDraw(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
145   //
146   // easy drawing of data, use '~' for abbreviation of '.fElements'
147   // example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
148  // sector: sector-number - only the specified sector will be drwawn
149   //         'A'/'C' or 'a'/'c' - side A/C will be drawn
150   //         'ALL' - whole TPC will be drawn, projected on one side
151   // cuts: specifies cuts
152   // drawOptions: draw options like 'same'
153   // writeDrawCommand: write the command, that is passed to TTree::Draw
154   //
155    TString drawStr(drawCommand);
156    TString sectorStr(sector);
157    sectorStr.ToUpper();
158    TString cutStr("");
159    TString drawOptionsStr("profcolz ");
160    TRandom rnd(0);
161    Int_t rndNumber = rnd.Integer(10000);
162    if (drawOptions && drawOptions != "")
163       drawOptionsStr += drawOptions;
164
165    if (sectorStr == "A") {
166       drawStr += ":gy.fElements:gx.fElements>>prof";
167       drawStr += rndNumber;
168       drawStr += "(330,-250,250,330,-250,250)";
169       cutStr += "(sector/18)%2==0 ";
170    }
171    else if  (sectorStr == "C") {
172       drawStr += ":gy.fElements:gx.fElements>>prof";
173       drawStr += rndNumber;
174       drawStr += "(330,-250,250,330,-250,250)";
175       cutStr += "(sector/18)%2==1 ";
176    }
177    else if  (sectorStr == "ALL") {
178       drawStr += ":gy.fElements:gx.fElements>>prof";
179       drawStr += rndNumber;
180       drawStr += "(330,-250,250,330,-250,250)";
181    }
182    else if (sectorStr.IsDigit()) {
183       Int_t isec = sectorStr.Atoi();
184       drawStr += ":rpad.fElements:row.fElements>>prof";
185       drawStr += rndNumber;
186       if (isec < 36 && isec >= 0)
187          drawStr += "(63,0,63,108,-54,54)";
188       else if (isec < 72 && isec >= 36)
189          drawStr += "(96,0,96,140,-70,70)";
190       else {
191          Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
192          return -1;
193       }
194       cutStr += "(sector==";
195       cutStr += isec;
196       cutStr += ") ";
197    }
198
199    if (cuts && cuts[0] != 0) {
200       if (cutStr.Length() != 0) cutStr += "&& ";
201       cutStr += "(";
202       cutStr += cuts;
203       cutStr += ")";
204    }
205    drawStr.ReplaceAll("~", ".fElements");
206    cutStr.ReplaceAll("~", ".fElements");
207    if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" <<  cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
208    return fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
209 }
210
211 Int_t AliTPCCalibViewer::EasyDraw(const char* drawCommand, Int_t sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
212   //
213   // easy drawing of data, use '~' for abbreviation of '.fElements'
214   // example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
215   // sector: sector-number - only the specified sector will be drwawn
216   // cuts: specifies cuts
217   // drawOptions: draw options like 'same'
218   // writeDrawCommand: write the command, that is passed to TTree::Draw
219   //
220    if (sector >= 0 && sector < 72) {
221       char sectorChr[3];
222       sprintf(sectorChr, "%i", sector);
223       return EasyDraw(drawCommand, sectorChr, cuts, drawOptions, writeDrawCommand);
224    }
225    Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
226    return -1;
227 }
228
229 //_____________________________________________________________________________
230 Int_t AliTPCCalibViewer::EasyDraw1D(const char* drawCommand, const char* sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
231   //
232   // easy drawing of data, use '~' for abbreviation of '.fElements'
233   // example: EasyDraw("CETmean~-CETmean_mean", "A", "(CETmean~-CETmean_mean)>0")
234   // sector: sector-number - the specified sector will be drwawn
235   //         'A'/'C' or 'a'/'c' - side A/C will be drawn
236   //         'ALL' - whole TPC will be drawn, projected on one side
237   // cuts: specifies cuts
238   // drawOptions: draw options like 'same'
239   // writeDrawCommand: write the command, that is passed to TTree::Draw
240   //
241
242    TString drawStr(drawCommand);
243    TString sectorStr(sector);
244    TString drawOptionsStr(drawOptions);
245    sectorStr.ToUpper();
246    TString cutStr("");
247
248    if (sectorStr == "A")
249       cutStr += "(sector/18)%2==0 ";
250    else if  (sectorStr == "C")
251       cutStr += "(sector/18)%2==1 ";
252    else if (sectorStr.IsDigit()) {
253       Int_t isec = sectorStr.Atoi();
254       if (isec < 0 || isec > 71) {
255          Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
256          return -1;
257       }
258       cutStr += "(sector==";
259       cutStr += isec;
260       cutStr += ") ";
261    }
262
263    if (cuts && cuts[0] != 0) {
264       if (cutStr.Length() != 0) cutStr += "&& ";
265       cutStr += "(";
266       cutStr += cuts;
267       cutStr += ")";
268    }
269
270    drawStr.ReplaceAll("~", ".fElements");
271    cutStr.ReplaceAll("~", ".fElements");
272    if (writeDrawCommand) std::cout << "fTree->Draw(\"" << drawStr << "\", \"" <<  cutStr << "\", \"" << drawOptionsStr << "\");" << std::endl;
273    return fTree->Draw(drawStr.Data(), cutStr.Data(), drawOptionsStr.Data());
274 }
275
276 Int_t AliTPCCalibViewer::EasyDraw1D(const char* drawCommand, Int_t sector, const char* cuts, const char* drawOptions, Bool_t writeDrawCommand) const {
277   //
278   // easy drawing of data, use '~' for abbreviation of '.fElements'
279   // example: EasyDraw("CETmean~-CETmean_mean", 34, "(CETmean~-CETmean_mean)>0")
280   // sector: sector-number - the specified sector will be drwawn
281   // cuts: specifies cuts
282   // drawOptions: draw options like 'same'
283   // writeDrawCommand: write the command, that is passed to TTree::Draw
284   //
285
286    if (sector >= 0 && sector < 72) {
287       char sectorChr[3];
288       sprintf(sectorChr, "%i", sector);
289       return EasyDraw1D(drawCommand, sectorChr, cuts, drawOptions, writeDrawCommand);
290    }
291   Error("EasyDraw","The TPC contains only sectors between 0 and 71.");
292   return -1;
293 }
294
295 //_____________________________________________________________________________
296 Int_t AliTPCCalibViewer::DrawHisto1D(const char* type, Int_t sector, TVectorF& nsigma, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM) const {
297   //
298   // draws a 1-dimensional histogram of 'type' for sector 'sector'
299   // TVectorF nsigma: Specifies, for which distances from the mean/median/LTM lines should be drawn, in units of sigma
300   // example: nsigma={2, 4, 6}: Three lines will be drawn, distance to mean/median/LTM: 2, 3 and 6 sigma
301   // plotMean, plotMedian, plotLTM: specifies, if mean, median and LTM should be drawn as lines into the histogram
302   //
303
304    TString typeStr(type);
305    TString sectorStr("sector==");
306    sectorStr += sector;
307
308    TCanvas* canvas = ((TCanvas*)gROOT->GetListOfCanvases()->Last());
309    Int_t oldOptStat = gStyle->GetOptStat();
310    gStyle->SetOptStat(0000000);
311
312    if (!canvas) {
313       canvas = new TCanvas();
314       fListOfObjectsToBeDeleted->Add(canvas);
315    }
316    
317    char c[500];
318    sprintf(c, "%s, sector: %i", type, sector);
319    TLegend * legend = new TLegend(.8,.6, .99, .99, c);
320    fListOfObjectsToBeDeleted->Add(legend);
321
322    Int_t nentries = fTree->Draw((typeStr+".fElements").Data(), sectorStr.Data(), "");
323    ((TH1F*)canvas->GetPrimitive("htemp"))->SetTitle("");
324       
325    //****************************************************************
326    //!!!!!!!!!!!!!!!!! Needs further investigaton !!!!!!!!!!!!!!!!!!!
327    //****************************************************************
328    //fListOfObjectsToBeDeleted->Add(canvas->GetPrimitive("htemp"));
329 /*
330   By default the temporary histogram created is called "htemp", but only in
331   the one dimensional Draw("e1") it contains the TTree's data points. For
332   a two dimensional Draw, the data is filled into a TGraph which is named
333   "Graph". They can be retrieved by calling
334     TH1F *htemp = (TH1F*)gPad->GetPrimitive("htemp"); // 1D
335     TGraph *graph = (TGraph*)gPad->GetPrimitive("Graph"); // 2D
336 */
337    
338    canvas->Update();
339    Double_t sigma = 0;
340
341    if (plotMean) {
342       fTree->Draw((typeStr+"_Mean").Data(), sectorStr.Data(), "goff");
343       Double_t lineX = fTree->GetV1()[0];
344       fTree->Draw((typeStr+"_RMS").Data(), sectorStr.Data(), "goff");
345       sigma = fTree->GetV1()[0];
346       TLine* line = new TLine(lineX, 0, lineX, canvas->GetUymax());
347       fListOfObjectsToBeDeleted->Add(line);
348       line->SetLineColor(kRed);
349       line->SetLineWidth(2);
350       line->SetLineStyle(1);
351       line->Draw();
352       sprintf(c, "Mean: %f", lineX);
353       legend->AddEntry(line, c, "l");
354
355       for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
356          TLine* linePlusSigma = new TLine(lineX+nsigma[i]*sigma, 0, lineX+nsigma[i]*sigma, canvas->GetUymax());
357          fListOfObjectsToBeDeleted->Add(linePlusSigma);
358          linePlusSigma->SetLineColor(kRed);
359          linePlusSigma->SetLineStyle(2+i);
360          linePlusSigma->Draw();
361    
362          TLine* lineMinusSigma = new TLine(lineX-nsigma[i]*sigma, 0, lineX-nsigma[i]*sigma, canvas->GetUymax());
363          fListOfObjectsToBeDeleted->Add(lineMinusSigma);
364          lineMinusSigma->SetLineColor(kRed);
365          lineMinusSigma->SetLineStyle(2+i);
366          lineMinusSigma->Draw();
367          sprintf(c, "%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i]*sigma));
368          std::cout << "nsigma-char*: " << c << std::endl;
369          legend->AddEntry(lineMinusSigma, c, "l");
370       }
371    }
372
373    if (plotMedian) {
374       fTree->Draw((typeStr+"_Median").Data(), sectorStr.Data(), "goff");
375       Double_t lineX = fTree->GetV1()[0];
376       fTree->Draw((typeStr+"_RMS").Data(), sectorStr.Data(), "goff");
377       sigma = fTree->GetV1()[0];
378       TLine* line = new TLine(lineX, 0, lineX, canvas->GetUymax());
379       fListOfObjectsToBeDeleted->Add(line);
380       line->SetLineColor(kBlue);
381       line->SetLineWidth(2);
382       line->SetLineStyle(1);
383       line->Draw();
384       sprintf(c, "Median: %f", lineX);
385       legend->AddEntry(line, c, "l");
386       
387       for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
388          TLine* linePlusSigma = new TLine(lineX+nsigma[i]*sigma, 0, lineX+nsigma[i]*sigma, canvas->GetUymax());
389          fListOfObjectsToBeDeleted->Add(linePlusSigma);
390          linePlusSigma->SetLineColor(kBlue);
391          linePlusSigma->SetLineStyle(2+i);
392          linePlusSigma->Draw();
393    
394          TLine* lineMinusSigma = new TLine(lineX-nsigma[i]*sigma, 0, lineX-nsigma[i]*sigma, canvas->GetUymax());
395          fListOfObjectsToBeDeleted->Add(lineMinusSigma);
396          lineMinusSigma->SetLineColor(kBlue);
397          lineMinusSigma->SetLineStyle(2+i);
398          lineMinusSigma->Draw();
399          sprintf(c, "%i #sigma = %f",(Int_t)(nsigma[i]), (Float_t)(nsigma[i]*sigma));
400          legend->AddEntry(lineMinusSigma, c, "l");
401       }
402    }
403    
404    if (plotLTM) {
405       fTree->Draw((typeStr+"_LTM").Data(), sectorStr.Data(), "goff");
406       Double_t lineX = fTree->GetV1()[0];
407       fTree->Draw((typeStr+"_RMS_LTM").Data(), sectorStr.Data(), "goff");
408       sigma = fTree->GetV1()[0];
409       TLine* line = new TLine(lineX, 0, lineX, canvas->GetUymax());
410       fListOfObjectsToBeDeleted->Add(line);
411       line->SetLineColor(kGreen+2);
412       line->SetLineWidth(2);
413       line->SetLineStyle(1);
414       line->Draw();
415       sprintf(c, "LTM: %f", lineX);
416       legend->AddEntry(line, c, "l");
417
418       for (Int_t i = 0; i < nsigma.GetNoElements(); i++) {
419          TLine* linePlusSigma = new TLine(lineX+nsigma[i]*sigma, 0, lineX+nsigma[i]*sigma, canvas->GetUymax());
420          fListOfObjectsToBeDeleted->Add(linePlusSigma);
421          linePlusSigma->SetLineColor(kGreen+2);
422          linePlusSigma->SetLineStyle(2+i);
423          linePlusSigma->Draw();
424    
425          TLine* lineMinusSigma = new TLine(lineX-nsigma[i]*sigma, 0, lineX-nsigma[i]*sigma, canvas->GetUymax());
426          fListOfObjectsToBeDeleted->Add(lineMinusSigma);
427          lineMinusSigma->SetLineColor(kGreen+2);
428          lineMinusSigma->SetLineStyle(2+i);
429          lineMinusSigma->Draw();
430          sprintf(c, "%i #sigma = %f", (Int_t)(nsigma[i]), (Float_t)(nsigma[i]*sigma));
431          legend->AddEntry(lineMinusSigma, c, "l");
432       }
433    }
434    legend->Draw();
435    gStyle->SetOptStat(oldOptStat);
436    return nentries;
437 }
438
439 //_____________________________________________________________________________
440 void AliTPCCalibViewer::SigmaCut(const char* type, Int_t sector, Float_t sigmaMax, Float_t sigmaStep, Bool_t plotMean, Bool_t plotMedian, Bool_t plotLTM) const {
441   //
442   // Creates a histogram, where you can see, how much of the data are inside sigma-intervals around the mean/median/LTM
443   // type: For which type of data the histogram is generated, e.g. 'CEQmean'
444   // sector: For which sector the histogram is generated
445   // sigmaMax: up to which sigma around the mean/median/LTM the histogram is generated
446   // sigmaStep: the binsize of the generated histogram
447   // plotMean/plotMedian/plotLTM: specifies where to put the center
448   //
449    Int_t oldOptStat = gStyle->GetOptStat();
450    gStyle->SetOptStat(0000000);
451
452    TString typeStr(type);
453    TString sectorStr("sector==");
454    sectorStr += sector;
455    Int_t entries = fTree->Draw((typeStr+".fElements").Data(), sectorStr.Data(), "goff");
456    char headline[500];
457    sprintf(headline, "%s in sector %i; Multiples of #sigma; Fraction of used pads", type, sector);
458    TH1F *histMean =   new TH1F("histMean",headline, (Int_t)(sigmaMax/sigmaStep+1), 0, sigmaMax+sigmaStep);
459    sprintf(headline, "%s in sector %i; Multiples of #sigma; Fraction of used pads", type, sector);
460    TH1F *histMedian = new TH1F("histMedian",headline, (Int_t)(sigmaMax/sigmaStep+1), 0, sigmaMax+sigmaStep);
461    sprintf(headline, "%s in sector %i; Multiples of #sigma; Fraction of used pads", type, sector);
462    TH1F *histLTM =    new TH1F("histLTM",headline, (Int_t)(sigmaMax/sigmaStep+1), 0, sigmaMax+sigmaStep);                                                                   
463    histMean->SetDirectory(0);
464    histMedian->SetDirectory(0);
465    histLTM->SetDirectory(0);
466    fListOfObjectsToBeDeleted->Add(histMean);
467    fListOfObjectsToBeDeleted->Add(histMedian);
468    fListOfObjectsToBeDeleted->Add(histLTM);
469
470    
471    // example-cut: sector==34 && TMath::Abs(CEQmean.fElements - CEQmean_Mean) < nsigma * CEQmean_RMS
472    for (Float_t nsigma = 0; nsigma <= sigmaMax; nsigma += sigmaStep) {
473      std::cout << "Calculating histograms,  step: " << (Int_t)(nsigma/sigmaStep) << " of: " << (Int_t)(sigmaMax/sigmaStep) << "\r" << std::flush;
474       char cuts[5000];
475       
476       if (plotMean) {
477          sprintf(cuts, "sector==%i && ( %s.fElements - %s_Median) < %f * %s_RMS", sector, type, type, nsigma, type );
478          sprintf(cuts,         "%s && (-%s.fElements + %s_Median) < %f * %s_RMS",   cuts, type, type, nsigma, type );
479          Float_t value = (Float_t)fTree->Draw((typeStr+".fElements").Data(), cuts, "goff")/entries;
480          histMean->Fill(nsigma, value);
481       }
482       if (plotMedian) {
483          sprintf(cuts, "sector==%i && ( %s.fElements - %s_Mean) < %f * %s_RMS", sector, type, type, nsigma, type );
484          sprintf(cuts,         "%s && (-%s.fElements + %s_Mean) < %f * %s_RMS",   cuts, type, type, nsigma, type );
485          Float_t value = (Float_t)fTree->Draw((typeStr+".fElements").Data(), cuts, "goff")/entries;
486          histMedian->Fill(nsigma, value);
487       }
488       if (plotLTM) {
489          sprintf(cuts, "sector==%i && ( %s.fElements - %s_LTM) < %f * %s_RMS_LTM", sector, type, type, nsigma, type );
490          sprintf(cuts,         "%s && (-%s.fElements + %s_LTM) < %f * %s_RMS_LTM",   cuts, type, type, nsigma, type );
491          Float_t value = (Float_t)fTree->Draw((typeStr+".fElements").Data(), cuts, "goff")/entries;
492          histLTM->Fill(nsigma, value);
493       }
494    }
495    
496    char c[500];
497    sprintf(c, "Sigma Cut");
498    TLegend * legend = new TLegend(.85,.8, .99, .99, c);
499    fListOfObjectsToBeDeleted->Add(legend);
500    
501    if (plotMean){
502       histMean->SetLineColor(kBlack);
503       sprintf(c, "Mean");
504       legend->AddEntry(histMean, c, "l");
505       histMean->Draw();
506    }
507    if (plotMedian){
508       histMedian->SetLineColor(kRed);
509       sprintf(c, "Median");
510       legend->AddEntry(histMedian, c, "l");
511       histMedian->Draw("same");
512    }
513    if (plotLTM){
514       histLTM->SetLineColor(kBlue);
515       sprintf(c, "LTM");
516       legend->AddEntry(histLTM, c, "l");
517       histLTM->Draw("same");
518    }   
519
520    legend->Draw();
521    gStyle->SetOptStat(oldOptStat);
522 }
523
524
525 //_____________________________________________________________________________
526 AliTPCCalPad* AliTPCCalibViewer::GetCalPad(const char* desiredData, char* cuts, char* calPadName) const {
527   //
528   // creates a AliTPCCalPad out of the 'desiredData'
529   // the functionality of EasyDraw1D is used
530   // calPadName specifies the name of the created AliTPCCalPad
531   //  - this takes a while -
532   //
533    TString drawStr(desiredData);
534    drawStr.Append(":channel~");
535    AliTPCCalPad * createdCalPad = new AliTPCCalPad(calPadName, calPadName);
536    Int_t entries = 0;
537    for (Int_t sec = 0; sec < 72; sec++) {
538       entries = EasyDraw1D(drawStr.Data(), (Int_t)sec, cuts, "goff");
539       for (Int_t i = 0; i < entries; i++) 
540          createdCalPad->GetCalROC(sec)->SetValue((UInt_t)(fTree->GetV2()[i]), (Float_t)(fTree->GetV1()[i]));
541    }
542    return createdCalPad;   
543 }
544
545 //_____________________________________________________________________________
546 AliTPCCalROC* AliTPCCalibViewer::GetCalROC(const char* desiredData, UInt_t sector, char* cuts) const {
547   //
548   // creates a AliTPCCalROC out of the desiredData
549   // the functionality of EasyDraw1D is used
550   // sector specifies the sector of the created AliTPCCalROC
551   //
552    TString drawStr(desiredData);
553    drawStr.Append(":channel~");
554    Int_t entries = EasyDraw1D(drawStr.Data(), (Int_t)sector, cuts, "goff");
555    AliTPCCalROC * createdROC = new AliTPCCalROC(sector);
556    for (Int_t i = 0; i < entries; i++) 
557       createdROC->SetValue((UInt_t)(fTree->GetV2()[i]), fTree->GetV1()[i]);
558    return createdROC;
559 }
560
561
562 TObjArray* AliTPCCalibViewer::GetListOfVariables(Bool_t printList) {
563   //
564   // scan the tree  - produces a list of available variables in the tree
565   // printList: print the list to the screen, after the scan is done
566   //
567    TObjArray* arr = new TObjArray();
568    TObjString* str = 0;
569    Int_t nentries = fTree->GetListOfBranches()->GetEntries();
570    for (Int_t i = 0; i < nentries; i++) {
571       str = new TObjString(fTree->GetListOfBranches()->At(i)->GetName());
572       str->String().ReplaceAll("_Median", "");
573       str->String().ReplaceAll("_Mean", "");
574       str->String().ReplaceAll("_RMS", "");
575       str->String().ReplaceAll("_LTM", "");
576       str->String().ReplaceAll("_OutlierCutted", "");
577       str->String().ReplaceAll(".", "");
578       if (!arr->FindObject(str) && 
579           !(str->String() == "channel" || str->String() == "gx" || str->String() == "gy" || 
580             str->String() == "lx" || str->String() == "ly" || str->String() == "pad" || 
581             str->String() == "row" || str->String() == "rpad" || str->String() == "sector"  ))
582          arr->Add(str);
583    }
584    arr->Sort();
585
586    if (printList) {
587       TIterator* iter = arr->MakeIterator();
588       iter->Reset();
589       TObjString* currentStr = 0;
590       while ( (currentStr = (TObjString*)(iter->Next())) ) {
591          std::cout << currentStr->GetString().Data() << std::endl;
592       }
593       delete iter;
594    }
595    return arr;
596 }
597
598
599 TObjArray* AliTPCCalibViewer::GetListOfNormalizationVariables(Bool_t printList) {
600   //
601   // produces a list of available variables for normalization in the tree
602   // printList: print the list to the screen, after the scan is done
603   //
604    TObjArray* arr = new TObjArray();
605    arr->Add(new TObjString("_Mean"));
606    arr->Add(new TObjString("_Mean_OutlierCutted"));
607    arr->Add(new TObjString("_Median"));
608    arr->Add(new TObjString("_Median_OutlierCutted"));
609    arr->Add(new TObjString("_LTM"));
610    arr->Add(new TObjString("_LTM_OutlierCutted"));
611    arr->Add(new TObjString("LFitIntern_4_8.fElements"));
612    arr->Add(new TObjString("GFitIntern_Lin.fElements"));
613    arr->Add(new TObjString("GFitIntern_Par.fElements"));
614
615    if (printList) {
616       TIterator* iter = arr->MakeIterator();
617       iter->Reset();
618       TObjString* currentStr = 0;
619       while ((currentStr = (TObjString*)(iter->Next()))) {
620          std::cout << currentStr->GetString().Data() << std::endl;
621       }
622       delete iter;
623    }
624    return arr;
625 }
626
627
628 TFriendElement* AliTPCCalibViewer::AddReferenceTree(const char* filename, const char* treename, const char* refname){
629   //
630   // add a reference tree to the current tree
631   // by default the treename is 'calPads' and the reference treename is 'R'
632   //
633    TFile *file = new TFile(filename);
634    fListOfObjectsToBeDeleted->Add(file);
635    TTree * tree = (TTree*)file->Get(treename);
636    return AddFriend(tree, refname);
637 }
638
639
640 TObjArray* AliTPCCalibViewer::GetArrayOfCalPads(){
641   //
642   // Returns a TObjArray with all AliTPCCalPads that are stored in the tree
643   //  - this takes a while - 
644   //
645    TObjArray *listOfCalPads = GetListOfVariables();
646    TObjArray *calPadsArray = new TObjArray();
647    Int_t numberOfCalPads = listOfCalPads->GetEntries();
648    for (Int_t i = 0; i < numberOfCalPads; i++) {
649      std::cout << "Creating calPad " << (i+1) << " of " << numberOfCalPads << "\r" << std::flush;
650       char* calPadName = (char*)((TObjString*)(listOfCalPads->At(i)))->GetString().Data();
651       TString drawCommand = ((TObjString*)(listOfCalPads->At(i)))->GetString();
652       drawCommand.Append("~");
653       AliTPCCalPad* calPad = GetCalPad(drawCommand.Data(), "", calPadName); 
654       calPadsArray->Add(calPad); 
655    }
656    std::cout << std::endl;
657    listOfCalPads->Delete();
658    delete listOfCalPads;
659    return calPadsArray;
660 }
661
662
663 void AliTPCCalibViewer::MakeTreeWithObjects(const char * fileName, TObjArray * array, const char * mapFileName) {
664   //
665   // Write tree with all available information
666   // im mapFileName is speciefied, the Map information are also written to the tree
667   // AliTPCCalPad-Objects are written directly to the tree, so that they can be accessd later on
668   // (does not work!!!)
669   //
670    AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
671
672    TObjArray* mapIROCs = 0;
673    TObjArray* mapOROCs = 0;
674    TVectorF *mapIROCArray = 0;
675    TVectorF *mapOROCArray = 0;
676    Int_t mapEntries = 0;
677    TString* mapNames = 0;
678    
679    if (mapFileName) {
680       TFile mapFile(mapFileName, "read");
681       
682       TList* listOfROCs = mapFile.GetListOfKeys();
683       mapEntries = listOfROCs->GetEntries()/2;
684       mapIROCs = new TObjArray(mapEntries*2);
685       mapOROCs = new TObjArray(mapEntries*2);
686       mapIROCArray = new TVectorF[mapEntries];
687       mapOROCArray = new TVectorF[mapEntries];
688       
689       mapNames = new TString[mapEntries];
690       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
691          TString ROCname(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
692          ROCname.Remove(ROCname.Length()-4, 4);
693          mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((ROCname + "IROC").Data()), ivalue);
694          mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((ROCname + "OROC").Data()), ivalue);
695          mapNames[ivalue].Append(ROCname);
696       }
697       
698       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
699          mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
700          mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
701       
702          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
703             (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
704          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
705             (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
706       }
707
708    } //  if (mapFileName)
709   
710    TTreeSRedirector cstream(fileName);
711    Int_t arrayEntries = array->GetEntries();
712    
713    // Read names of AliTPCCalPads and save them in names[]
714    TString* names = new TString[arrayEntries];
715    for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
716       names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
717
718    for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
719       
720       TVectorF *vectorArray = new TVectorF[arrayEntries];
721       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
722          vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
723             
724       
725       //
726       // fill vectors of variable per pad
727       //
728       TVectorF *posArray = new TVectorF[8];
729       for (Int_t ivalue = 0; ivalue < 8; ivalue++)
730          posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
731
732       Float_t posG[3] = {0};
733       Float_t posL[3] = {0};
734       Int_t ichannel = 0;
735       for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
736          for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
737             tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
738             tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
739             posArray[0][ichannel] = irow;
740             posArray[1][ichannel] = ipad;
741             posArray[2][ichannel] = posL[0];
742             posArray[3][ichannel] = posL[1];
743             posArray[4][ichannel] = posG[0];
744             posArray[5][ichannel] = posG[1];
745             posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
746             posArray[7][ichannel] = ichannel;
747             
748             // loop over array containing AliTPCCalPads
749             for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
750                AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
751                AliTPCCalROC* calROC = calPad->GetCalROC(isector);
752                if (calROC)
753                   (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
754                else
755                   (vectorArray[ivalue])[ichannel] = 0;
756             }
757             ichannel++;
758          }
759       }
760       AliTPCCalROC dummyROC(0);
761       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
762          AliTPCCalROC *roc = ((AliTPCCalPad*)array->At(ivalue))->GetCalROC(isector);
763          if (!roc) roc = &dummyROC;
764          cstream << "calPads" <<
765             (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
766          cstream << "calPads" << 
767             (Char_t*)((names[ivalue] + "Pad.=").Data()) << roc;
768       }
769
770       if (mapFileName) {
771          for  (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
772             if (isector < 36)
773                cstream << "calPads" <<
774                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
775             else
776                cstream << "calPads" <<
777                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
778          }
779       }
780       
781       cstream << "calPads" <<
782          "sector=" << isector;
783
784       cstream << "calPads" <<
785          "row.=" << &posArray[0] <<
786          "pad.=" << &posArray[1] <<
787          "lx.=" << &posArray[2] <<
788          "ly.=" << &posArray[3] <<
789          "gx.=" << &posArray[4] <<
790          "gy.=" << &posArray[5] <<
791          "rpad.=" << &posArray[6] <<
792          "channel.=" << &posArray[7];
793
794       cstream << "calPads" <<
795          "\n";
796
797       delete[] posArray;
798       delete[] vectorArray;
799    } //for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++)
800
801    delete[] names;
802    if (mapFileName) {
803       delete mapIROCs;
804       delete mapOROCs;
805       delete[] mapIROCArray;
806       delete[] mapOROCArray;
807       delete[] mapNames;
808    }
809 }
810
811 void AliTPCCalibViewer::MakeTree(const char * fileName, TObjArray * array, const char * mapFileName, AliTPCCalPad* outlierPad, Float_t ltmFraction) {
812   //
813   // Write a tree with all available information
814   // im mapFileName is speciefied, the Map information are also written to the tree
815   // pads specified in outlierPad are not used for calculating statistics
816   //  - the same function as AliTPCCalPad::MakeTree - 
817   //
818    AliTPCROC* tpcROCinstance = AliTPCROC::Instance();
819
820    TObjArray* mapIROCs = 0;
821    TObjArray* mapOROCs = 0;
822    TVectorF *mapIROCArray = 0;
823    TVectorF *mapOROCArray = 0;
824    Int_t mapEntries = 0;
825    TString* mapNames = 0;
826    
827    if (mapFileName) {
828       TFile mapFile(mapFileName, "read");
829       
830       TList* listOfROCs = mapFile.GetListOfKeys();
831       mapEntries = listOfROCs->GetEntries()/2;
832       mapIROCs = new TObjArray(mapEntries*2);
833       mapOROCs = new TObjArray(mapEntries*2);
834       mapIROCArray = new TVectorF[mapEntries];
835       mapOROCArray = new TVectorF[mapEntries];
836       
837       mapNames = new TString[mapEntries];
838       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
839          TString ROCname(((TKey*)(listOfROCs->At(ivalue*2)))->GetName());
840          ROCname.Remove(ROCname.Length()-4, 4);
841          mapIROCs->AddAt((AliTPCCalROC*)mapFile.Get((ROCname + "IROC").Data()), ivalue);
842          mapOROCs->AddAt((AliTPCCalROC*)mapFile.Get((ROCname + "OROC").Data()), ivalue);
843          mapNames[ivalue].Append(ROCname);
844       }
845       
846       for (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
847          mapIROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(0));
848          mapOROCArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(36));
849       
850          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(0); ichannel++)
851             (mapIROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapIROCs->At(ivalue)))->GetValue(ichannel);
852          for (UInt_t ichannel = 0; ichannel < tpcROCinstance->GetNChannels(36); ichannel++)
853             (mapOROCArray[ivalue])[ichannel] = ((AliTPCCalROC*)(mapOROCs->At(ivalue)))->GetValue(ichannel);
854       }
855
856    } //  if (mapFileName)
857   
858    TTreeSRedirector cstream(fileName);
859    Int_t arrayEntries = array->GetEntries();
860    
861    TString* names = new TString[arrayEntries];
862    for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
863       names[ivalue].Append(((AliTPCCalPad*)array->At(ivalue))->GetName());
864
865    for (UInt_t isector = 0; isector < tpcROCinstance->GetNSectors(); isector++) {
866       //
867       // get statistic for given sector
868       //
869       TVectorF median(arrayEntries);
870       TVectorF mean(arrayEntries);
871       TVectorF rms(arrayEntries);
872       TVectorF ltm(arrayEntries);
873       TVectorF ltmrms(arrayEntries);
874       TVectorF medianWithOut(arrayEntries);
875       TVectorF meanWithOut(arrayEntries);
876       TVectorF rmsWithOut(arrayEntries);
877       TVectorF ltmWithOut(arrayEntries);
878       TVectorF ltmrmsWithOut(arrayEntries);
879       
880       TVectorF *vectorArray = new TVectorF[arrayEntries];
881       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++)
882          vectorArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
883       
884       for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
885          AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
886          AliTPCCalROC* calROC = calPad->GetCalROC(isector);
887          AliTPCCalROC* outlierROC = 0;
888          if (outlierPad) outlierROC = outlierPad->GetCalROC(isector);
889          if (calROC) {
890             median[ivalue] = calROC->GetMedian();
891             mean[ivalue] = calROC->GetMean();
892             rms[ivalue] = calROC->GetRMS();
893             Double_t ltmrmsValue = 0;
894             ltm[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction);
895             ltmrms[ivalue] = ltmrmsValue;
896             if (outlierROC) {
897                medianWithOut[ivalue] = calROC->GetMedian(outlierROC);
898                meanWithOut[ivalue] = calROC->GetMean(outlierROC);
899                rmsWithOut[ivalue] = calROC->GetRMS(outlierROC);
900                ltmrmsValue = 0;
901                ltmWithOut[ivalue] = calROC->GetLTM(&ltmrmsValue, ltmFraction, outlierROC);
902                ltmrmsWithOut[ivalue] = ltmrmsValue;
903             }
904          }
905          else {
906             median[ivalue] = 0.;
907             mean[ivalue] = 0.;
908             rms[ivalue] = 0.;
909             ltm[ivalue] = 0.;
910             ltmrms[ivalue] = 0.;
911             medianWithOut[ivalue] = 0.;
912             meanWithOut[ivalue] = 0.;
913             rmsWithOut[ivalue] = 0.;
914             ltmWithOut[ivalue] = 0.;
915             ltmrmsWithOut[ivalue] = 0.;
916          }
917       }
918       
919       //
920       // fill vectors of variable per pad
921       //
922       TVectorF *posArray = new TVectorF[8];
923       for (Int_t ivalue = 0; ivalue < 8; ivalue++)
924          posArray[ivalue].ResizeTo(tpcROCinstance->GetNChannels(isector));
925
926       Float_t posG[3] = {0};
927       Float_t posL[3] = {0};
928       Int_t ichannel = 0;
929       for (UInt_t irow = 0; irow < tpcROCinstance->GetNRows(isector); irow++) {
930          for (UInt_t ipad = 0; ipad < tpcROCinstance->GetNPads(isector, irow); ipad++) {
931             tpcROCinstance->GetPositionLocal(isector, irow, ipad, posL);
932             tpcROCinstance->GetPositionGlobal(isector, irow, ipad, posG);
933             posArray[0][ichannel] = irow;
934             posArray[1][ichannel] = ipad;
935             posArray[2][ichannel] = posL[0];
936             posArray[3][ichannel] = posL[1];
937             posArray[4][ichannel] = posG[0];
938             posArray[5][ichannel] = posG[1];
939             posArray[6][ichannel] = (Int_t)(ipad - (Double_t)(tpcROCinstance->GetNPads(isector, irow))/2);
940             posArray[7][ichannel] = ichannel;
941             
942             // loop over array containing AliTPCCalPads
943             for (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
944                AliTPCCalPad* calPad = (AliTPCCalPad*) array->At(ivalue);
945                AliTPCCalROC* calROC = calPad->GetCalROC(isector);
946                if (calROC)
947                   (vectorArray[ivalue])[ichannel] = calROC->GetValue(irow, ipad);
948                else
949                   (vectorArray[ivalue])[ichannel] = 0;
950             }
951             ichannel++;
952          }
953       }
954       
955       cstream << "calPads" <<
956          "sector=" << isector;
957       
958       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
959          cstream << "calPads" <<
960             (Char_t*)((names[ivalue] + "_Median=").Data()) << median[ivalue] <<
961             (Char_t*)((names[ivalue] + "_Mean=").Data()) << mean[ivalue] <<
962             (Char_t*)((names[ivalue] + "_RMS=").Data()) << rms[ivalue] <<
963             (Char_t*)((names[ivalue] + "_LTM=").Data()) << ltm[ivalue] <<
964             (Char_t*)((names[ivalue] + "_RMS_LTM=").Data()) << ltmrms[ivalue];
965          if (outlierPad) {
966             cstream << "calPads" <<
967                (Char_t*)((names[ivalue] + "_Median_OutlierCutted=").Data()) << medianWithOut[ivalue] <<
968                (Char_t*)((names[ivalue] + "_Mean_OutlierCutted=").Data()) << meanWithOut[ivalue] <<
969                (Char_t*)((names[ivalue] + "_RMS_OutlierCutted=").Data()) << rmsWithOut[ivalue] <<
970                (Char_t*)((names[ivalue] + "_LTM_OutlierCutted=").Data()) << ltmWithOut[ivalue] <<
971                (Char_t*)((names[ivalue] + "_RMS_LTM_OutlierCutted=").Data()) << ltmrmsWithOut[ivalue];
972          }
973       }
974
975       for  (Int_t ivalue = 0; ivalue < arrayEntries; ivalue++) {
976          cstream << "calPads" <<
977             (Char_t*)((names[ivalue] + ".=").Data()) << &vectorArray[ivalue];
978       }
979
980       if (mapFileName) {
981          for  (Int_t ivalue = 0; ivalue < mapEntries; ivalue++) {
982             if (isector < 36)
983                cstream << "calPads" <<
984                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapIROCArray[ivalue];
985             else
986                cstream << "calPads" <<
987                   (Char_t*)((mapNames[ivalue] + ".=").Data()) << &mapOROCArray[ivalue];
988          }
989       }
990
991       cstream << "calPads" <<
992          "row.=" << &posArray[0] <<
993          "pad.=" << &posArray[1] <<
994          "lx.=" << &posArray[2] <<
995          "ly.=" << &posArray[3] <<
996          "gx.=" << &posArray[4] <<
997          "gy.=" << &posArray[5] <<
998          "rpad.=" << &posArray[6] <<
999          "channel.=" << &posArray[7];
1000          
1001       cstream << "calPads" <<
1002          "\n";
1003
1004       delete[] posArray;
1005       delete[] vectorArray;
1006    }
1007    
1008
1009    delete[] names;
1010    if (mapFileName) {
1011       delete mapIROCs;
1012       delete mapOROCs;
1013       delete[] mapIROCArray;
1014       delete[] mapOROCArray;
1015       delete[] mapNames;
1016    }
1017 }
1018