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