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