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