]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDQADataMakerRec.cxx
Removing shadowed variables
[u/mrichter/AliRoot.git] / TRD / AliTRDQADataMakerRec.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 /* $Id$ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 //                                                                        //
20 //  Produces the data needed to calculate the quality assurance.          //
21 //  All data must be mergeable objects.                                   //
22 //                                                                        //
23 //  Author:                                                               //
24 //    Sylwester Radomski (radomski@physi.uni-heidelberg.de)               //
25 //                                                                        //
26 ////////////////////////////////////////////////////////////////////////////
27
28 // --- ROOT system ---
29 #include <TClonesArray.h>
30 #include <TFile.h> 
31 #include <TH1D.h> 
32 #include <TH2D.h>
33 #include <TH3D.h>
34 #include <TProfile.h>
35 #include <TF1.h>
36 #include <TCanvas.h>
37
38 // --- AliRoot header files ---
39 #include "AliESDEvent.h"
40 #include "AliLog.h"
41 #include "AliTRDcluster.h"
42 #include "AliTRDQADataMakerRec.h"
43 #include "AliTRDgeometry.h"
44 #include "AliTRDdataArrayI.h"
45 #include "AliTRDrawStreamTB.h"
46
47 #include "AliQAChecker.h"
48
49 ClassImp(AliTRDQADataMakerRec)
50
51 //____________________________________________________________________________ 
52   AliTRDQADataMakerRec::AliTRDQADataMakerRec() : 
53   AliQADataMakerRec(AliQA::GetDetName(AliQA::kTRD), "TRD Quality Assurance Data Maker")
54 {
55   //
56   // Default constructor
57 }
58
59 //____________________________________________________________________________ 
60 AliTRDQADataMakerRec::AliTRDQADataMakerRec(const AliTRDQADataMakerRec& qadm) :
61   AliQADataMakerRec()
62 {
63   //
64   // Copy constructor 
65   //
66
67   SetName((const char*)qadm.GetName()) ; 
68   SetTitle((const char*)qadm.GetTitle()); 
69
70 }
71
72 //__________________________________________________________________
73 AliTRDQADataMakerRec& AliTRDQADataMakerRec::operator=(const AliTRDQADataMakerRec& qadm)
74 {
75   //
76   // Equal operator.
77   //
78
79   this->~AliTRDQADataMakerRec();
80   new(this) AliTRDQADataMakerRec(qadm);
81   return *this;
82
83 }
84
85 //____________________________________________________________________________ 
86 void AliTRDQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray * list)
87 {
88   //
89   // Detector specific actions at end of cycle
90   //
91   //TStopwatch watch;
92   //watch.Start();
93   
94   AliInfo("End of TRD cycle");
95   
96   if (task == AliQA::kRECPOINTS) {
97     
98     TH1D *hist = new TH1D("fitHist", "", 200, -0.5, 199.5);
99     //list->Print();
100     
101     // fill detector map;
102     for(int i=0; i<540; i++) {
103       Double_t v = ((TH1D*)list->At(0))->GetBinContent(i+1);
104       Int_t sm = i/30;
105       Int_t det = i%30;
106
107       TH2D *detMap = (TH2D*)list->At(87);
108       Int_t bin = detMap->FindBin(sm, det);
109       detMap->SetBinContent(bin, v);
110     }
111
112
113     // Rec points full chambers
114     for (Int_t i=0; i<540; i++) {
115       
116       //AliInfo(Form("I = %d", i));
117
118       //TH1D *h = ((TH2D*)list->At(1))->ProjectionY(Form("qaTRD_recPoints_amp_%d",i), i+1, i+1);
119       hist->Reset();
120       for(Int_t b=1; b<hist->GetXaxis()->GetNbins()-1; b++) {
121         Double_t xvalue = hist->GetBinCenter(b);
122         Int_t bin = ((TH2D*)list->At(1))->FindBin(i,xvalue);
123         Double_t value =  ((TH2D*)list->At(1))->GetBinContent(bin);
124         hist->SetBinContent(b, value);
125       }
126       
127       //printf("Sum = %d %f\n", i, hist->GetSum());
128       if (hist->GetSum() < 100) continue; // chamber not present
129       
130       hist->Fit("landau", "q0", "goff", 10, 180);
131       TF1 *fit = hist->GetFunction("landau");
132       ((TH1D*)list->At(12))->Fill(fit->GetParameter(1));
133       ((TH1D*)list->At(13))->Fill(fit->GetParameter(2));
134     }
135     
136     // time-bin by time-bin sm by sm
137     for(Int_t i=0; i<18; i++) { // loop over super-modules
138       
139       for(Int_t j=0; j<35; j++) { // loop over time bins
140         
141         //TH1D *h =  ((TH3D*)list->At(10))->ProjectionZ(Form("ampTime_%d",i), i+1, i+1, j+1, j+1);     
142         hist->Reset();
143         for(Int_t b=1; b<hist->GetXaxis()->GetNbins()-1; b++) {
144           Double_t xvalue = hist->GetBinCenter(b);
145           Double_t svalue = 0;
146           
147           for(Int_t det=i*30; det<(i+1)*30; det++) { // loop over detectors
148             Int_t bin = ((TH3D*)list->At(10))->FindBin(det,j,xvalue);
149             Double_t value =  ((TH3D*)list->At(10))->GetBinContent(bin);
150             svalue += value;
151           }
152           //printf("v = %f\n", value);
153           hist->SetBinContent(b, svalue);
154         }
155         
156         if (hist->GetSum() < 100) continue;
157         //printf("fitting %d %d %f\n", i, j, hist->GetSum());
158         
159         hist->Fit("landau", "q0", "goff", 10, 180);
160         TF1 *fit = hist->GetFunction("landau");
161         
162         TH1D *h1 = (TH1D*)list->At(14+18+i);
163         Int_t bin = h1->FindBin(j);
164         // printf("%d %d %d\n", det, j, bin);
165         h1->SetBinContent(bin, TMath::Abs(fit->GetParameter(1)));
166       }
167     }
168       
169
170     // time-bin by time-bin chamber by chamber
171     for (Int_t i=0; i<540; i++) {
172       
173       //TH1D *test = ((TH3D*)list->At(10))->ProjectionZ(Form("ampTime_%d",i), i+1, i+1, 0, 35);     
174       //if (test->GetSum() < 100) continue;
175       
176       //AliInfo(Form("fitting det = %d", i));
177       
178       for(Int_t j=0; j<35; j++) {
179         
180         //TH1D *h =  ((TH3D*)list->At(10))->ProjectionZ(Form("ampTime_%d",i), i+1, i+1, j+1, j+1);     
181         hist->Reset();
182         for(Int_t b=1; b<hist->GetXaxis()->GetNbins()-1; b++) {
183           Double_t xvalue = hist->GetBinCenter(b);
184           Int_t bin = ((TH3D*)list->At(10))->FindBin(i,j,xvalue);
185           Double_t value =  ((TH3D*)list->At(10))->GetBinContent(bin);
186           //printf("v = %f\n", value);
187           hist->SetBinContent(b, value);
188         }
189         
190         if (hist->GetSum() < 100) continue;
191         //printf("fitting %d %d %f\n", i, j, hist->GetSum());
192         
193         hist->Fit("landau", "q0", "goff", 10, 180);
194         TF1 *fit = hist->GetFunction("landau");
195         
196         Int_t sm = i/30;
197         Int_t det = i%30;
198         TH2D *h2 = (TH2D*)list->At(14+sm);
199         Int_t bin = h2->FindBin(det,j);
200         // printf("%d %d %d\n", det, j, bin);
201         h2->SetBinContent(bin, TMath::Abs(fit->GetParameter(1)));
202         h2->SetBinError(bin,fit->GetParError(1));
203       }
204     }
205     
206     if (hist) delete hist;
207   }
208   
209   //////////////////////////
210   // const Int_t knbits = 6;
211   // const char *suf[knbits] = {"TPCi", "TPCo", "TPCz", "TRDo", "TRDr", "TRDz"};
212   //const char *sufRatio[4] = {"TRDrTRDo", "TRDoTPCo", "TRDrTPCo", "TRDzTPCo"};
213
214   if (task == AliQA::kESDS) {
215     
216     const Int_t knRatio = 4;
217     const Int_t kN[knRatio] = {4,3,4,5};
218     const Int_t kD[knRatio] = {3,1,1,3}; 
219     
220     // create ratios
221     for(Int_t type=0; type<2; type++) {
222       for(Int_t i=0; i<knRatio; i++) {
223
224         TH1D *ratio = (TH1D*)list->At(19 + 2*i + type);
225         TH1D *histN = (TH1D*)list->At(3 + 2*kN[i] + type);
226         TH1D *histD = (TH1D*)list->At(3 + 2*kD[i] + type);
227
228         BuildRatio(ratio, histN, histD);
229         //ratio->Reset();
230         //ratio->Add(histN);
231         //ratio->Divide(histD);
232       }
233     }
234
235     // ratio for the fraction of electrons per stack
236     TH1D *histN = (TH1D*)list->At(33);
237     TH1D *histD = (TH1D*)list->At(32);
238     TH1D *ratio = (TH1D*)list->At(34);
239     BuildRatio(ratio, histN, histD);
240   }
241
242   
243   // call the checker
244   AliQAChecker::Instance()->Run(AliQA::kTRD, task, list) ;    
245 }
246
247 //____________________________________________________________________________ 
248 void AliTRDQADataMakerRec::InitESDs()
249 {
250   //
251   // Create ESDs histograms in ESDs subdir
252   //
253
254   const Int_t kNhist = 36+5+4;
255
256   TH1 *hist[kNhist];
257   Int_t histoCounter = -1 ;
258
259   hist[++histoCounter] = new TH1D("qaTRD_esd_ntracks", ";Number of tracks", 300, -0.5, 299.5);
260   hist[++histoCounter] = new TH1D("qaTRD_esd_sector", ";Sector", 18, -0.5, 17.7);
261   hist[++histoCounter] = new TH1D("qaTRD_esd_bits", ";Bits", 64, -0.5, 63.5);
262
263   const Int_t knbits = 6;
264   const char *suf[knbits] = {"TPCi", "TPCo", "TPCz", "TRDo", "TRDr", "TRDz"};
265
266   // histo = 3
267   for(Int_t i=0; i<knbits; i++) { 
268     hist[++histoCounter] = new TH1D(Form("qaTRD_esd_pt%s",suf[i]), ";p_{T} (GeV/c);", 100, 0, 10);
269     hist[++histoCounter] = new TH1D(Form("qaTRD_esd_trdz%s", suf[i]), ";z (cm)", 200, -400, 400); 
270   }
271
272   hist[++histoCounter] = new TH1D("qaTRD_esd_clsTRDo", "TRDo;number of clusters", 180, -0.5, 179.5);;
273   hist[++histoCounter] = new TH1D("qaTRD_esd_clsTRDr", "TRDr;number of clusters", 180, -0.5, 179.5);;
274   hist[++histoCounter] = new TH1D("qaTRD_esd_clsTRDz", "TRDz;number of clusters", 180, -0.5, 179.5);;
275   //hist[++histoCounter] = new TH1D("qaTRD_esd_clsRatio", ";cluster ratio", 100, 0., 1.3);;
276
277   hist[++histoCounter] = new TH2D("qaTRD_esd_sigMom", ";momentum (GeV/c);signal", 100, 0, 5, 200, 0, 1e3);
278
279   // end of cycle plots (hist 19)
280   const char *sufRatio[4] = {"TRDrTRDo", "TRDoTPCo", "TRDrTPCo", "TRDzTPCo"};
281
282   for(int i=0; i<4; i++) {
283     hist[++histoCounter] = new TH1D(Form("qaTRD_esd_pt%s",sufRatio[i]), 
284                                     Form("Efficiency in Pt %s;p_{T};eff", sufRatio[i]),
285                                     100, 0, 10);
286
287     hist[++histoCounter] = new TH1D(Form("qaTRD_esd_trdz%s",sufRatio[i]), 
288                                     Form("Efficiency in Z %s;z (cm);eff", sufRatio[i]),
289                                     200, -400, 400);
290   }
291
292   // 27 - 31
293   hist[27] = new TH1D("qaTRD_esd_quality", ";quality", 120, 0, 12);
294   hist[28] = new TH1D("qaTRD_esd_budget", ";NN", 110, -1000, 100);
295   hist[29] = new TH1D("qaTRD_esd_chi2", ";chi2", 200, 0, 100);
296   hist[30] = new TH1D("qaTRD_esd_timeBin", ";time bin", 7, -0.5, 6.5);
297   hist[31] = new TH1D("qaTRD_esd_pidQuality", "pid Quality", 7, -0.5, 6.5);
298
299   // stack by stack electron identyfication
300   hist[32] = new TH1D("qaTRD_esd_tracksStack", "number of all tracks;stack", 90, -0.5, 89.5);
301   hist[33] = new TH1D("qaTRD_esd_electronStack", "number of electron tracks;stack", 90, -0.5, 89.5);
302   hist[34] = new TH1D("qaTRD_esd_elRatioStack", "fraction of electron tracks;stack", 90, -0.5, 89.5);
303   hist[35] = new TH1D("qaTRD_esd_thetaOut", ";tan(theta);", 100, -1, 1);
304   
305   const char *partType[5] = {"Electron", "Muon", "Pion", "Kaon", "Proton"}; 
306
307   for(Int_t i=0; i<AliPID::kSPECIES; i++)
308     hist[36+i] = new TH1D(Form("qaTRD_esd_pid%d",i),
309                           Form("%s;probability",partType[i]), 100, 0, 1);
310  
311   // dE/dX vs momentum in three regions
312   const char *zoneName[4] = {"total charge", "ampilification range", "plateau", "TR range"};
313  
314   // prepare the scale from 0.1 to 10 GeV
315   const Int_t nscalex= 50;
316   Double_t scalex[nscalex+1];
317   Double_t dd = (TMath::Log(10) - TMath::Log(0.5)) / nscalex;
318   for(Int_t ix=0; ix<nscalex+1; ix++) {
319     scalex[ix] = 0.5 * TMath::Exp(dd * ix);
320   }
321
322   const Int_t nscaley = 50;
323   Double_t scaley[nscaley+1];
324   for(Int_t iy=0; iy<nscaley+1; iy++) {
325     scaley[iy] = iy * (3e3/nscaley);
326   }
327     
328
329   for(Int_t i=0; i<4; i++) {
330     hist[41+i] = new TH2D(Form("qaTRD_esd_signalPzone_%d",i), 
331                           Form("%s;momentum (GeV/c);singal (a.u.)", zoneName[i]),
332                           nscalex, scalex, nscaley, scaley);
333   }
334
335   for(Int_t i=0; i<kNhist; i++) {
336     //hist[i]->Sumw2();
337     Add2ESDsList(hist[i], i);
338   }
339
340 }
341
342 //____________________________________________________________________________ 
343 void AliTRDQADataMakerRec::InitRecPoints()
344 {
345   //
346   // Create Reconstructed Points histograms in RecPoints subdir
347   //
348
349   const Int_t kNhist = 14 + 4 * 18 + 2;
350   TH1 *hist[kNhist];
351
352   hist[0] = new TH1D("qaTRD_recPoints_det", ";Detector ID of the cluster", 540, -0.5, 539.5);
353   hist[1] = new TH2D("qaTRD_recPoints_amp", ";Amplitude", 540, -0.5, 539, 200, -0.5, 199.5);
354   hist[2] = new TH1D("qaTRD_recPoints_npad", ";Number of Pads", 12, -0.5, 11.5);
355
356   hist[3] = new TH1D("qaTRD_recPoints_dist2", ";residuals [2pad]", 100, -1, 1);
357   hist[4] = new TH1D("qaTRD_recPoints_dist3", ";residuals [3pad]", 100, -1, 1);
358   hist[5] = new TH1D("qaTRD_recPoints_dist4", ";residuals [4pad]", 100, -1, 1);
359   hist[6] = new TH1D("qaTRD_recPoints_dist5", ";residuals [5pad]", 100, -1, 1);
360
361   hist[7] = new TH2D("qaTRD_recPoints_rowCol", ";row;col", 16, -0.5, 15.5, 145, -0.5, 144.5);
362   hist[8] = new TH1D("qaTRD_recPoints_time", ";time bin", 35, -0.5, 34.5);
363   hist[9] = new TH1D("qaTRD_recPoints_nCls", ";number of clusters", 500, -0.5, 499.5);
364
365   hist[10] = new TH3D("qaTRD_recPoints_sigTime", ";chamber;time bin;signal", 
366                       540, -0.5, 539.5, 35, -0.5, 34.5, 200, -0.5, 199.5);
367   hist[11] = new TProfile("qaTRD_recPoints_prf", ";distance;center of gravity"
368                          , 120, -0.6, 0.6, -1.2, 1.2, "");
369
370   hist[12] = new TH1D("qaTRD_recPoints_ampMPV", ";amplitude MPV", 200, 0, 200);
371   hist[13] = new TH1D("qaTRD_recPoints_ampSigma", ";amplitude Sigma", 200, 0, 200); 
372   
373   // chamber by chamber
374   for(Int_t i=0; i<18; i++) {
375     hist[14+i] = new TH2D(Form("qaTRD_recPoints_sigTime_sm%d",i), Form("sm%d;det;time bin"), 
376                         30, -0.5, 29.5, 35, -0.5, 34.5);
377     hist[14+i]->SetMinimum(20);
378     hist[14+i]->SetMaximum(40);
379   }
380  
381   // time bin by time bin sm-by-sm
382   for(Int_t i=0; i<18; i++) {
383     hist[14+18+i] = new TH1D(Form("qaTRD_recPoints_sigTimeShape_sm%d", i), 
384                              Form("sm%d;time bin;signal"),
385                              35, -0.5, 34.5);
386
387     hist[14+18+i]->SetMaximum(120);    
388   }
389
390   // str = 50
391   for(Int_t i=0; i<18; i++) {
392     hist[50+i] = new TH1D(Form("qaTRD_recPoints_nCls_sm%d",i),
393                           Form("sm%d;time bin;number of clusters",i),
394                           35, -0.5, 34.5);
395   }
396
397   // str = 68
398   for(Int_t i=0; i<18; i++) {
399     hist[68+i] = new TH1D(Form("qaTRD_recPoints_totalCharge_sm%d", i),
400                           Form("sm%d;time bin;total charge", i),
401                           35, -0.5, 34.5);
402   }
403
404   hist[86] = new TH1D("qaTRD_recPoints_signal", ";amplitude", 200, -0.5, 199.5);
405   hist[87] = new TH2D("qaTRD_recPoints_detMap", ";sm;chamber", 18, -0.5, 17.5, 30, -0.5, 29.5);
406
407
408   for(Int_t i=0; i<kNhist; i++) {
409     //hist[i]->Sumw2();
410     Add2RecPointsList(hist[i], i);
411   }
412 }
413
414 //____________________________________________________________________________ 
415 void AliTRDQADataMakerRec::InitRaws()
416 {
417   //
418   // create Raws histograms in Raws subdir
419   //
420
421   const Int_t kSM = 18;
422   //const Int_t kNCh = 540;
423   const Int_t kNhist = 4+kSM;
424   TH1D *hist[kNhist];
425
426   // four histograms to be published
427   hist[0] = new TH1D("qaTRD_raws_det", ";detector", 540, -0.5, 539.5);
428   hist[1] = new TH1D("qaTRD_raws_sig", ";signal", 100, -0.5, 99.5);
429   hist[2] = new TH1D("qaTRD_raws_timeBin", ";time bin", 40, -0.5, 39.5); 
430   hist[3] = new TH1D("qaTRD_raws_smId", ";supermodule", 18, -0.5, 17.5);
431   //
432   
433   // one double per MCM (not published)
434   const Int_t kNMCM = 30 * 8 * 16;
435   for(Int_t i=0; i<kSM; i++)
436     hist[4+i] = new TH1D(Form("qaTRD_raws_sm%d",i),"",kNMCM, -0.5, kNMCM-0.5); 
437   
438   // register
439   for(Int_t i=0; i<kNhist; i++) {
440     //hist[i]->Sumw2();
441     Add2RawsList(hist[i], i);
442   }
443
444 }
445
446 //____________________________________________________________________________
447 void AliTRDQADataMakerRec::MakeESDs(AliESDEvent * esd)
448 {
449   //
450   // Make QA data from ESDs
451   //
452
453   Int_t nTracks = esd->GetNumberOfTracks();
454   GetESDsData(0)->Fill(nTracks);
455
456   // track loop
457   for (Int_t i=0; i<nTracks; i++) {
458
459     AliESDtrack *track = esd->GetTrack(i);
460     const AliExternalTrackParam *paramOut = track->GetOuterParam();
461     const AliExternalTrackParam *paramIn = track->GetInnerParam();
462
463     // long track ..
464     if (!paramIn) continue;
465     if (!paramOut) continue;
466
467     // not a kink
468     if (track->GetKinkIndex(0) > 0) continue; 
469
470     Double_t extZ = GetExtZ(paramIn);
471     if (TMath::Abs(extZ) > 320) continue; // acceptance cut
472
473     // .. in the acceptance
474     Int_t sector = GetSector(paramOut->GetAlpha());
475     Int_t stack = GetStack(paramOut);
476
477     UInt_t u = 1;
478     UInt_t status = track->GetStatus();
479     for(Int_t bit=0; bit<32; bit++) 
480       if (u<<bit & status) GetESDsData(2)->Fill(bit);
481
482     const Int_t knbits = 6; 
483     Int_t bit[6] = {0,0,0,0,0,0};    
484     bit[0] = status & AliESDtrack::kTPCin;
485     bit[1] = status & AliESDtrack::kTPCout;
486     bit[2] = (status & AliESDtrack::kTPCout) && !(status & AliESDtrack::kTRDout);
487     bit[3] = status & AliESDtrack::kTRDout;
488     bit[4] = status & AliESDtrack::kTRDrefit;
489     bit[5] = (status & AliESDtrack::kTRDout) && !(status & AliESDtrack::kTRDrefit);
490
491     // transverse momentum
492     //const Double_t *val = paramOut->GetParameter(); // parameters at the Outer plane
493     Double_t pt = paramOut->Pt(); //1./TMath::Abs(val[4]);
494
495     for(Int_t b=0; b<knbits; b++) {
496       if (bit[b]) {
497         GetESDsData(2*b+3)->Fill(pt); 
498         GetESDsData(2*b+4)->Fill(extZ);
499       }
500     }
501
502     // clusters
503     for(Int_t b=0; b<3; b++) 
504       if (bit[3+b]) GetESDsData(b+15)->Fill(track->GetTRDncls0());
505
506     // refitted only
507     if (!bit[4]) continue;
508
509     //fQuality->Fill(track->GetTRDQuality());
510     //fBudget->Fill(track->GetTRDBudget());
511     //fSignal->Fill(track->GetTRDsignal());
512
513     GetESDsData(1)->Fill(sector);
514     GetESDsData(18)->Fill(track->GetP(), track->GetTRDsignal());
515
516     GetESDsData(27)->Fill(track->GetTRDQuality());
517     GetESDsData(28)->Fill(track->GetTRDBudget());
518     GetESDsData(29)->Fill(track->GetTRDchi2());
519     GetESDsData(30)->Fill(track->GetTRDTimBin(0));
520     GetESDsData(31)->Fill(track->GetTRDpidQuality());
521     
522     
523     // dedx
524     for(Int_t k=0; k<4; ++k) {
525       Double_t dedx = 0;
526       for(Int_t j=0; j<6; j++) {
527         dedx += track->GetTRDslice(j, k-1);
528       }
529       GetESDsData(41+i)->Fill(paramOut->GetP(), dedx/6.);
530     }
531
532     // probabilities
533     if (status & AliESDtrack::kTRDpid) {
534       for(Int_t k=0; k<AliPID::kSPECIES; ++k) 
535         GetESDsData(36+k)->Fill(track->GetTRDpid(k));
536     }
537
538     // probabilities uniformity
539     if (track->GetTRDpidQuality() < 6) continue;
540     GetESDsData(35)->Fill(paramOut->GetZ()/paramOut->GetX());
541     
542     Int_t idx = 5 * sector + stack;
543     GetESDsData(32)->Fill(idx); // all tracks
544     if (track->GetTRDpid(AliPID::kElectron) > 0.9) 
545       GetESDsData(33)->Fill(idx); // electrons only
546
547     
548
549     /*
550     hist[27] = new TH1D("qaTRD_esd_quality", ";quality", 120, 0, 12);
551     hist[28] = new TH1D("qaTRD_esd_budget", ";NN", 110, -1000, 100);
552     hist[29] = new TH1D("qaTRD_esd_chi2", ";chi2", 300, 0, 100);
553     hist[30] = new TH1D("qaTRD_esd_timeBin", 7, -0.5, 6.5);
554     hist[31] = new TH1D("qaTRD_esd_pidQuality", 7, -0.5, 6.5);
555     */
556
557     /*
558     // PID only
559     if (status & AliESDtrack::kTRDpid) {
560
561       for(Int_t l=0; l<6; l++) fTime->Fill(track->GetTRDTimBin(l));
562
563       // fill pid histograms
564       Double_t trdr0 = 0; //, tpcr0 = 0;
565       Int_t trdBestPid = 5; //, tpcBestPid = 5;  // charged
566       const Double_t kminPidValue = 0.9;
567
568       //Double_t pp[5];
569       //track->GetTPCpid(pp); // ESD inconsequence
570
571       for(Int_t pid=0; pid<5; pid++) {
572         
573         trdr0 += track->GetTRDpid(pid);
574         //tpcr0 += pp[pid];
575         
576         fTrdPID[pid]->Fill(track->GetTRDpid(pid));
577         //fTpcPID[pid]->Fill(pp[pid]);
578         
579         if (track->GetTRDpid(pid) > kminPidValue) trdBestPid = pid;
580         //if (pp[pid] > kminPidValue) tpcBestPid = pid;
581       }
582
583       fTrdPID[5]->Fill(trdr0); // check unitarity
584       fTrdSigMomPID[trdBestPid]->Fill(track->GetP(), track->GetTRDsignal());
585
586       //fTpcPID[5]->Fill(tpcr0); // check unitarity
587       //fTpcSigMomPID[tpcBestPid]->Fill(track->GetP(), track->GetTPCsignal());
588     }
589     */
590
591   }
592
593 }
594
595 //______________________________________________________________________________
596 Int_t AliTRDQADataMakerRec::GetSector(Double_t alpha) const 
597 {
598   //
599   // Gets the sector number 
600   //
601
602   Double_t size = TMath::DegToRad() * 20.; // shall use TRDgeo
603   if (alpha < 0) alpha += 2*TMath::Pi();
604   Int_t sector = (Int_t)(alpha/size);
605   return sector;
606
607 }
608 //______________________________________________________________________________
609
610 Int_t AliTRDQADataMakerRec::GetStack(const AliExternalTrackParam *paramOut) const
611 {
612   //
613   // calculates the stack the track is in
614   //
615   
616   const Double_t L = -0.9;
617   const Double_t W = (2*L)/5;
618
619   Double_t tan = paramOut->GetZ() / paramOut->GetX();
620   Double_t pos = (tan - L) / W;
621   return (Int_t) pos;
622 }
623
624 //______________________________________________________________________________
625 Double_t AliTRDQADataMakerRec::GetExtZ(const AliExternalTrackParam *in) const 
626 {
627   //
628   // Returns the Z position at the entry to TRD
629   // using parameters from the TPC in
630   //
631
632   const Double_t kX0 = 300;
633
634   Double_t x = in->GetX();
635   const Double_t *par = in->GetParameter();
636   Double_t theta = par[3];
637   Double_t z = in->GetZ();
638
639   Double_t zz = z + (kX0-x) * TMath::Tan(theta);
640   return zz;
641
642 }
643
644 //____________________________________________________________________________
645 void AliTRDQADataMakerRec::MakeRaws(AliRawReader* rawReader)
646 {
647   //
648   // Makes QA data from raw data
649   //
650
651   // 157
652   // T9 -- T10
653
654   //const Int_t kSM = 18;
655   //const Int_t kROC = 30;
656   const Int_t kROB = 8;
657   //const Int_t kLayer = 6;
658   //const Int_t kStack = 5;
659   const Int_t kMCM = 16;
660   //  const Int_t kADC = 22;
661
662   //AliTRDrawStreamBase::SetRawStreamVersion("TB");
663   AliTRDrawStreamBase *raw = AliTRDrawStreamBase::GetRawStream(rawReader);
664   AliDebug(2,Form("Stream version: %s", raw->IsA()->GetName()));
665
666   while (raw->Next()) {
667
668     GetRawsData(0)->Fill(raw->GetDet());
669
670     // possibly needs changes with the new reader !!
671     Int_t *sig = raw->GetSignals();
672     for(Int_t i=0; i<3; i++) GetRawsData(1)->Fill(sig[i]);
673     // ---
674
675     GetRawsData(2)->Fill(raw->GetTimeBin());
676
677     // calculate the index;
678     Int_t sm = raw->GetSM();
679     Int_t roc = raw->GetROC();
680     Int_t rob = raw->GetROB();
681     Int_t mcm = raw->GetMCM();
682     //Int_t adc = raw->GetADC();
683
684     //Int_t index = roc * (kROB*kMCM*kADC) + rob * (kMCM*kADC) + mcm * kADC + adc;
685     Int_t  index = roc * (kROB*kMCM) + rob * kMCM + mcm;
686     GetRawsData(3)->Fill(sm);
687     GetRawsData(4+sm)->Fill(index);
688   }
689
690   delete raw;
691 }
692
693 //____________________________________________________________________________
694 void AliTRDQADataMakerRec::MakeRecPoints(TTree * clustersTree)
695 {
696   //  
697   // Makes data from RecPoints
698   // 
699   
700   //  Info("MakeRecPoints", "making");
701
702   Int_t nsize = Int_t(clustersTree->GetTotBytes() / (sizeof(AliTRDcluster))); 
703   TObjArray *clusterArray = new TObjArray(nsize+1000); 
704
705   TBranch *branch = clustersTree->GetBranch("TRDcluster");
706   if (!branch) {
707     AliError("Can't get the branch !");
708     return;
709   }
710   branch->SetAddress(&clusterArray); 
711
712   // Loop through all entries in the tree
713   Int_t nEntries   = (Int_t) clustersTree->GetEntries();
714   Int_t nbytes     = 0;
715   AliTRDcluster *c = 0;
716   Int_t nDet[540];
717   for (Int_t i=0; i<540; i++) nDet[i] = 0;
718
719   for (Int_t iEntry = 0; iEntry < nEntries; iEntry++) {    
720
721     // Import the tree
722     nbytes += clustersTree->GetEvent(iEntry);  
723
724     // Get the number of points in the detector
725     Int_t nCluster = clusterArray->GetEntriesFast();  
726
727     // Loop through all TRD digits
728     for (Int_t iCluster = 0; iCluster < nCluster; iCluster++) { 
729       c = (AliTRDcluster *) clusterArray->UncheckedAt(iCluster);
730
731       Int_t iDet = c->GetDetector();
732       nDet[iDet]++;
733       GetRecPointsData(0)->Fill(iDet);
734       GetRecPointsData(86)->Fill(c->GetQ());
735       GetRecPointsData(1)->Fill(iDet, c->GetQ());
736       GetRecPointsData(2)->Fill(c->GetNPads());
737       if (c->GetNPads() < 6)
738         GetRecPointsData(1+c->GetNPads())->Fill(c->GetCenter());
739
740       //if (c->GetPadTime() < 5)
741       ((TH2D*)GetRecPointsData(7))->Fill(c->GetPadRow(), c->GetPadCol());
742       GetRecPointsData(8)->Fill(c->GetPadTime());
743
744       ((TH3D*)GetRecPointsData(10))->Fill(iDet, c->GetPadTime(), c->GetQ());
745       
746       Int_t iSM = iDet / 30;
747       GetRecPointsData(50+iSM)->Fill(c->GetPadTime());
748       GetRecPointsData(68+iSM)->Fill(c->GetPadTime(), c->GetQ());
749
750       // PRF for 2pad
751       //if (c->GetNPads() == 2) {
752       Short_t *sig = c->GetSignals();
753       Double_t frac = -10;
754
755       if (sig[0] == 0 && sig[1] == 0 && sig[2] == 0 && sig[5] == 0 && sig[6] == 0) 
756         frac = 1. * sig[4] / (sig[3] + sig[4]);
757
758       if (sig[0] == 0 && sig[1] == 0 && sig[4] == 0 && sig[5] == 0 && sig[6] == 0)
759         frac = -1. * sig[2] / (sig[2] + sig[3]);
760
761       if (frac > -10)  ((TProfile*)GetRecPointsData(11))->Fill(c->GetCenter(), frac);
762         
763       //}
764     }
765   }
766
767   for(Int_t i=0; i<540; i++) 
768     if (nDet[i] > 0) GetRecPointsData(9)->Fill(nDet[i]);
769
770   delete clusterArray;
771
772 }
773
774 //____________________________________________________________________________ 
775 void AliTRDQADataMakerRec::StartOfDetectorCycle()
776 {
777   //
778   // Detector specific actions at start of cycle
779   //
780
781 }
782
783 //__________________________________________________________________________
784 Int_t AliTRDQADataMakerRec::CheckPointer(TObject *obj, const char *name) 
785 {
786   //
787   // Checks initialization of pointers
788   //
789
790   if (!obj) AliWarning(Form("null pointer: %s", name));
791   return !!obj;
792 }
793 //__________________________________________________________________________
794 void AliTRDQADataMakerRec::BuildRatio(TH1D *ratio, TH1D *histN, TH1D*histD) {
795   //
796   // Calculate the ratio of two histograms 
797   // error are calculated assuming the histos have the same counts
798   //
799
800   // calclate
801
802   Int_t nbins = histN->GetXaxis()->GetNbins();
803   for(Int_t i=1; i<nbins+2; i++) {
804     
805     Double_t valueN = histN->GetBinContent(i);
806     Double_t valueD = histD->GetBinContent(i);
807     
808     if (valueD < 1) {
809       ratio->SetBinContent(i, 0);
810       ratio->SetBinError(i, 0);
811       continue;
812     }
813
814     Double_t eps = (valueN < valueD-valueN)? valueN : valueD-valueN;
815     
816     ratio->SetBinContent(i, valueN/valueD);
817     ratio->SetBinError(i, TMath::Sqrt(eps)/valueD);
818   }
819
820   // style
821   ratio->SetMinimum(-0.1);
822   ratio->SetMaximum(1.1);
823   ratio->SetMarkerStyle(20);
824 }
825 //__________________________________________________________________________