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