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