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