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