]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQADataMakerRec.cxx
3db65f9f454228660daf0c60a0577794c4e8fcf6
[u/mrichter/AliRoot.git] / ITS / AliITSQADataMakerRec.cxx
1 /**************************************************************************
2  * Copyright(c) 2007-2009, 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 //  Checks the quality assurance 
20 //  by comparing with reference data
21 //  contained in a DB
22 //  -------------------------------------------------------------
23 //  W. Ferrarese + P. Cerello Feb 2008
24 //  INFN Torino
25 //  Melinda Siciliano Aug 2008
26
27
28 // --- ROOT system ---
29 #include <TH2.h>
30 #include <TTree.h>
31 // --- Standard library ---
32
33 // --- AliRoot header files ---
34 #include "AliITSQADataMakerRec.h"
35 #include "AliITSQASPDDataMakerRec.h"
36 #include "AliITSQASDDDataMakerRec.h"
37 #include "AliITSQASSDDataMakerRec.h"
38 #include "AliLog.h"
39 #include "AliQAv1.h"
40 #include "AliQAChecker.h"
41 #include "AliITSQAChecker.h"
42 #include "AliITSRecPoint.h"
43 #include "AliRawReader.h"
44 #include "AliESDEvent.h"
45 #include "AliESDtrack.h"
46 #include "AliESDVertex.h"
47 #include "AliMultiplicity.h"
48 #include "AliITSgeomTGeo.h"
49
50 ClassImp(AliITSQADataMakerRec)
51
52 //____________________________________________________________________________ 
53 AliITSQADataMakerRec::AliITSQADataMakerRec(Bool_t kMode, Short_t subDet, Short_t ldc) :
54 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kITS), "ITS Quality Assurance Data Maker"),
55 fkOnline(kMode),
56 fSubDetector(subDet),
57 fLDC(ldc),
58 fRunNumber(0),
59 fSPDDataMaker(NULL),
60 fSDDDataMaker(NULL),
61 fSSDDataMaker(NULL)
62
63 {
64   //ctor used to discriminate OnLine-Offline analysis
65   if(fSubDetector < 0 || fSubDetector > 3) {
66         AliError("Error: fSubDetector number out of range; return\n");
67   }
68
69   // Initialization for RAW data 
70   if(fSubDetector == 0 || fSubDetector == 1) {
71     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Create SPD DataMakerRec\n");
72         fSPDDataMaker = new AliITSQASPDDataMakerRec(this,fkOnline);
73   }
74   if(fSubDetector == 0 || fSubDetector == 2) {
75     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Create SDD DataMakerRec\n");
76         fSDDDataMaker = new AliITSQASDDDataMakerRec(this,fkOnline);
77   }
78   if(fSubDetector == 0 || fSubDetector == 3) {
79     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Create SSD DataMakerRec\n");
80         fSSDDataMaker = new AliITSQASSDDataMakerRec(this,fkOnline);
81   }
82 }
83
84 //____________________________________________________________________________ 
85 AliITSQADataMakerRec::~AliITSQADataMakerRec(){
86   // destructor
87   if(fSPDDataMaker)delete fSPDDataMaker;
88   if(fSDDDataMaker)delete fSDDDataMaker;
89   if(fSSDDataMaker)delete fSSDDataMaker;
90 }
91
92 //____________________________________________________________________________ 
93 AliITSQADataMakerRec::AliITSQADataMakerRec(const AliITSQADataMakerRec& qadm) :
94 AliQADataMakerRec(),
95 fkOnline(qadm.fkOnline),
96 fSubDetector(qadm.fSubDetector),
97 fLDC(qadm.fLDC),
98 fRunNumber(0),
99 fSPDDataMaker(NULL),
100 fSDDDataMaker(NULL),
101 fSSDDataMaker(NULL)
102
103 {
104   //copy ctor 
105   SetName((const char*)qadm.GetName()) ; 
106   SetTitle((const char*)qadm.GetTitle());
107 }
108
109 //__________________________________________________________________
110 AliITSQADataMakerRec& AliITSQADataMakerRec::operator = (const AliITSQADataMakerRec& qac )
111 {
112   // Equal operator.
113   this->~AliITSQADataMakerRec();
114   new(this) AliITSQADataMakerRec(qac);
115   return *this;
116 }
117
118 //____________________________________________________________________________ 
119 void AliITSQADataMakerRec::StartOfDetectorCycle()
120 {
121   //Detector specific actions at start of cycle
122   AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM::Start of ITS Cycle\n");
123   if(fSubDetector == 0 || fSubDetector == 1) fSPDDataMaker->StartOfDetectorCycle();
124   if(fSubDetector == 0 || fSubDetector == 2) fSDDDataMaker->StartOfDetectorCycle();
125   if(fSubDetector == 0 || fSubDetector == 3) fSSDDataMaker->StartOfDetectorCycle();
126 }
127
128 //____________________________________________________________________________ 
129 void AliITSQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray** list)
130 {
131   // launch the QA checking
132
133   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
134     //SetEventSpecie(specie) ; 
135     //AliQAv1::Instance()->SetEventSpecie(specie);
136     if(AliQAv1::Instance()->IsEventSpecieSet(specie)){
137     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list[specie])\n"); 
138     if(fSubDetector == 0 || fSubDetector == 1) fSPDDataMaker->EndOfDetectorCycle(task, list[/*GetEventSpecie()*/specie]);
139     if(fSubDetector == 0 || fSubDetector == 2) fSDDDataMaker->EndOfDetectorCycle(task, list[/*GetEventSpecie()*/specie]);
140     if(fSubDetector == 0 || fSubDetector == 3) fSSDDataMaker->EndOfDetectorCycle(task, list[/*GetEventSpecie()*/specie]);
141   
142   
143     AliQAChecker *qac = AliQAChecker::Instance();
144     AliITSQAChecker *qacb = (AliITSQAChecker *) qac->GetDetQAChecker(0);
145     Int_t subdet=GetSubDet();
146     qacb->SetSubDet(subdet);
147  
148     if(subdet== 0 ){
149       qacb->SetTaskOffset(fSPDDataMaker->GetOffset(task), fSDDDataMaker->GetOffset(task), fSSDDataMaker->GetOffset(task)); //Setting the offset for the QAChecker list
150     }
151     else
152       if(subdet!=0){
153         Int_t offset=GetDetTaskOffset(subdet, task);
154         qacb->SetDetTaskOffset(subdet,offset);
155       }
156
157     qac->Run( AliQAv1::kITS , task, list);
158  
159     }
160     
161   }
162 }
163
164 //____________________________________________________________________________ 
165 void AliITSQADataMakerRec::EndOfDetectorCycle(const char * /*fgDataName*/)
166 {
167   //eventually used for different  AliQAChecker::Instance()->Run
168 }
169
170 //____________________________________________________________________________ 
171 void AliITSQADataMakerRec::InitRaws()
172 {  
173
174   //if(fRawsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries()) return;
175         
176
177   if(fSubDetector == 0 || fSubDetector == 1) {
178     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SPD InitRaws\n");
179     fSPDDataMaker->InitRaws();
180   }
181   if(fSubDetector == 0 || fSubDetector == 2) {
182     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SDD InitRaws\n");
183     
184     fSDDDataMaker->SetOffset(AliQAv1::kRAWS, fRawsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries(),AliRecoParam::AConvert(fEventSpecie));
185     fSDDDataMaker->InitRaws();
186   }
187   if(fSubDetector == 0 || fSubDetector == 3) {
188     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SSD InitRaws\n");
189     
190     fSSDDataMaker->SetOffset(AliQAv1::kRAWS, fRawsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries(),AliRecoParam::AConvert(fEventSpecie));
191     fSSDDataMaker->InitRaws();
192   }
193 }
194
195 //____________________________________________________________________________
196 void AliITSQADataMakerRec::MakeRaws(AliRawReader* rawReader)
197
198   // Fill QA for RAW   
199   //return ; 
200
201   SetRunNumber(rawReader->GetRunNumber());
202
203   if(fSubDetector == 0 || fSubDetector == 1)  {
204     fSPDDataMaker->MakeRaws(rawReader) ; 
205   }
206   
207   if(fSubDetector == 0 || fSubDetector == 2) {
208     fSDDDataMaker->MakeRaws(rawReader) ; 
209   }
210
211   if(fSubDetector == 0 || fSubDetector == 3) fSSDDataMaker->MakeRaws(rawReader);
212
213 }
214
215 //____________________________________________________________________________ 
216 void AliITSQADataMakerRec::InitDigits()
217 {
218   
219   // Initialization for DIGITS
220   if(fSubDetector == 0 || fSubDetector == 1) {
221     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SPD InitDigitss\n");
222
223     fSPDDataMaker->InitDigits();
224   }
225   if(fSubDetector == 0 || fSubDetector == 2) {
226     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SDD InitDigits\n");
227     fSDDDataMaker->SetOffset(AliQAv1::kDIGITSR, fDigitsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries(),AliRecoParam::AConvert(fEventSpecie));
228
229     fSDDDataMaker->InitDigits();
230   }
231   if(fSubDetector == 0 || fSubDetector == 3) {
232     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SSD InitDigits\n");
233     fSSDDataMaker->SetOffset(AliQAv1::kDIGITSR, fDigitsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries(),AliRecoParam::AConvert(fEventSpecie));
234
235     fSSDDataMaker->InitDigits();
236   }
237 }
238
239 //____________________________________________________________________________ 
240 void AliITSQADataMakerRec::MakeDigits(TTree * digitsTree)
241 {
242
243   
244   // Fill QA for recpoints
245   if(fSubDetector == 0 || fSubDetector == 1) {
246     fSPDDataMaker->MakeDigits(digitsTree) ; 
247   }
248   
249   if(fSubDetector == 0 || fSubDetector == 2) {
250     fSDDDataMaker->MakeDigits(digitsTree) ; 
251     
252   }
253   
254   if(fSubDetector == 0 || fSubDetector == 3) fSSDDataMaker->MakeDigits(digitsTree);
255 }
256
257 //____________________________________________________________________________ 
258 void AliITSQADataMakerRec::InitRecPoints()
259 {
260
261   // Initialization for RECPOINTS
262
263
264   //if(fRecPointsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries()) return;
265         if(fSubDetector == 0 || fSubDetector == 1) {
266                 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SPD InitRecPoints\n");
267                 fSPDDataMaker->InitRecPoints();
268         }
269         if(fSubDetector == 0 || fSubDetector == 2) {
270                 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SDD InitRecPoints\n");
271                 fSDDDataMaker->SetOffset(AliQAv1::kRECPOINTS, fRecPointsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries(), AliRecoParam::AConvert(fEventSpecie));
272                 fSDDDataMaker->InitRecPoints();
273         }
274         if(fSubDetector == 0 || fSubDetector == 3) {
275                 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SSD InitRecPoints\n");
276                 fSSDDataMaker->SetOffset(AliQAv1::kRECPOINTS, fRecPointsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries(),AliRecoParam::AConvert(fEventSpecie));
277                 fSSDDataMaker->InitRecPoints();
278         }
279
280
281         if(fSubDetector == 0){
282           Int_t offset = fRecPointsQAList [AliRecoParam::AConvert(fEventSpecie)]->GetEntries();
283           const Bool_t expert   = kTRUE ; 
284           const Bool_t image    = kTRUE ; 
285           Char_t name[50];
286           Char_t title[50];
287           TH2F**hPhiEta = new TH2F*[6];
288           for (Int_t iLay=0;iLay<6;iLay++) {
289             sprintf(name,"Phi_vs_Eta_ITS_Layer%d",iLay+1);
290             sprintf(title,"Phi vs Eta - ITS Layer %d",iLay+1);
291             hPhiEta[iLay]=new TH2F(name,title,30,-1.5,1.5,200,0.,2*TMath::Pi());
292             hPhiEta[iLay]->GetXaxis()->SetTitle("Pseudorapidity");
293             hPhiEta[iLay]->GetYaxis()->SetTitle("#varphi [rad]");
294             Add2RecPointsList((new TH2F(*hPhiEta[iLay])), iLay + offset, !expert, image);
295             
296             delete hPhiEta[iLay];
297           }
298           
299         }
300         
301 }
302
303 //____________________________________________________________________________ 
304 void AliITSQADataMakerRec::MakeRecPoints(TTree * clustersTree)
305 {
306  
307   // Fill QA for recpoints
308   if(fSubDetector == 0 || fSubDetector == 1) {
309     fSPDDataMaker->MakeRecPoints(clustersTree) ; 
310   }
311     
312   if(fSubDetector == 0 || fSubDetector == 2) {
313     fSDDDataMaker->MakeRecPoints(clustersTree) ; 
314   }
315   
316   if(fSubDetector == 0 || fSubDetector == 3) fSSDDataMaker->MakeRecPoints(clustersTree);
317
318
319
320   if(fSubDetector == 0){
321     // Check id histograms already created for this Event Specie
322     TBranch *branchRecP = clustersTree->GetBranch("ITSRecPoints");
323     if (!branchRecP) {
324       AliError("can't get the branch with the ITS clusters !");
325       return;
326     }
327     
328     Int_t offset = fRecPointsQAList [AliRecoParam::AConvert(fEventSpecie)]->GetEntries();
329     Float_t cluGlo[3] = {0.,0.,0.};
330     Int_t lay, lad, det; 
331     static TClonesArray statRecpoints("AliITSRecPoint") ;
332     TClonesArray *recpoints = &statRecpoints;
333     branchRecP->SetAddress(&recpoints);
334     // Fill QA for recpoints
335     for(Int_t module=0; module<clustersTree->GetEntries();module++){
336       //  AliInfo(Form("Module %d\n",module));
337       branchRecP->GetEvent(module);
338       AliITSgeomTGeo::GetModuleId(module, lay, lad, det);
339       for(Int_t j=0;j<recpoints->GetEntries();j++){
340         AliITSRecPoint *rcp = (AliITSRecPoint*)recpoints->At(j);    
341         //Check id histograms already created for this Event Specie
342         rcp->GetGlobalXYZ(cluGlo);
343         Double_t rad=TMath::Sqrt(cluGlo[0]*cluGlo[0]+cluGlo[1]*cluGlo[1]+cluGlo[2]*cluGlo[2]);
344         Double_t phi= TMath::Pi() + TMath::ATan2(-cluGlo[1],-cluGlo[0]);
345         Double_t theta = TMath::ACos(cluGlo[2]/rad);
346         Double_t eta = 100.;
347         if(AreEqual(rad,0.) == kFALSE) {
348           if(theta<=1.e-14){ eta=30.; }
349           else { eta = -TMath::Log(TMath::Tan(theta/2.));}
350         }
351         //      printf("=========================>hlt   rcp->GetLayer() = %d \n",rcp->GetLayer());
352         (GetRecPointsData( rcp->GetLayer() + offset - 6))->Fill(eta,phi);
353       }
354     }
355   }
356   
357 }
358
359 //____________________________________________________________________________ 
360 void AliITSQADataMakerRec::FillRecPoint(AliITSRecPoint rcp)
361 {
362
363         // Fill QA for recpoints
364         Float_t cluGlo[3] = {0.,0.,0.};
365         Int_t offset = fRecPointsQAList [AliRecoParam::AConvert(fEventSpecie)]->GetEntries();
366   // Check id histograms already created for this Event Specie
367         rcp.GetGlobalXYZ(cluGlo);
368         Double_t rad=TMath::Sqrt(cluGlo[0]*cluGlo[0]+cluGlo[1]*cluGlo[1]+cluGlo[2]*cluGlo[2]);
369         Double_t phi= TMath::Pi() + TMath::ATan2(-cluGlo[1],-cluGlo[0]);
370         Double_t theta = TMath::ACos(cluGlo[2]/rad);
371         Double_t eta = 100.;
372         if(AreEqual(rad,0.)==kFALSE) {
373           if(theta<=1.e-14){eta=30.;}
374           else    {eta = -TMath::Log(TMath::Tan(theta/2.));}
375         }
376         (GetRecPointsData( rcp.GetLayer() + offset - 6))->Fill(eta,phi);        
377
378 }
379
380 //____________________________________________________________________________ 
381 TH2F *AliITSQADataMakerRec::GetITSGlobalHisto(Int_t layer)
382 {
383
384         Int_t offset = fRecPointsQAList [AliRecoParam::AConvert(fEventSpecie)]->GetEntries();
385         return ((TH2F *) GetRecPointsData( layer + offset - 6));//local distribution
386 }
387
388 //____________________________________________________________________________ 
389 void AliITSQADataMakerRec::InitESDs()
390 {
391
392   // Create ESDs histograms in ESDs subdir
393
394   Bool_t expertHistogram = kTRUE;
395
396   TH1F *hESDClustersMI = 
397     new TH1F("hESDClustersMI", "N ITS clusters per track (MI); N clusters; Counts",
398              7, -0.5, 6.5);
399   hESDClustersMI->Sumw2();
400   hESDClustersMI->SetMinimum(0);
401   Add2ESDsList(hESDClustersMI, 0);
402
403   TH1F *hESDClusterMapMI =
404     new TH1F("hESDClusterMapMI", "N tracks with point Layer (MI); Layer; N tracks",
405              6, -0.5, 5.5);
406   hESDClusterMapMI->Sumw2();
407   hESDClusterMapMI->SetMinimum(0);
408   Add2ESDsList(hESDClusterMapMI, 1, expertHistogram);
409
410   TH1F *hESDClustersSA = 
411     new TH1F("hESDClustersSA", "N ITS clusters per track (SA); N clusters; Counts",
412              7, -0.5, 6.5);
413   hESDClustersSA->Sumw2();
414   hESDClustersSA->SetMinimum(0);
415   Add2ESDsList(hESDClustersSA, 2);
416
417   TH1F *hESDClusterMapSA =
418     new TH1F("hESDClusterMapSA", "N tracks with point Layer (SA); Layer; N tracks",
419              6, -0.5, 5.5);
420   hESDClusterMapSA->Sumw2();
421   hESDClusterMapSA->SetMinimum(0);
422   Add2ESDsList(hESDClusterMapSA, 3, expertHistogram);
423
424   TH1F *hSPDVertexX = 
425     new TH1F("hSPDVertexX","SPD Vertex x; x [cm]; N events",
426              10000,-2,2);
427   hSPDVertexX->Sumw2();
428   Add2ESDsList(hSPDVertexX, 4);
429
430   TH1F *hSPDVertexY = 
431     new TH1F("hSPDVertexY","SPD Vertex y; y [cm]; N events",
432              10000,-2,2);
433   hSPDVertexY->Sumw2();
434   Add2ESDsList(hSPDVertexY, 5);
435
436   TH1F *hSPDVertexZ = 
437     new TH1F("hSPDVertexZ","SPD Vertex Z; z [cm]; N events",
438              10000,-20,20);
439   hSPDVertexZ->Sumw2();
440   Add2ESDsList(hSPDVertexZ, 6);
441
442   TH1F *hSPDVertexContrOverMult =
443     new TH1F("hSPDVertexContrOverMult","SPD Vertex: contributors / multiplicity; N contributors / SPD multiplicity; N events",
444              100,-4,20);
445   hSPDVertexContrOverMult->Sumw2();
446   Add2ESDsList(hSPDVertexContrOverMult, 7, expertHistogram);
447
448   TH1F *hTrkVertexX = 
449     new TH1F("hTrkVertexX","ITS+TPC Trk Vertex x; x [cm]; N events",
450              10000,-2,2);
451   hTrkVertexX->Sumw2();
452   Add2ESDsList(hTrkVertexX, 8, expertHistogram);
453
454   TH1F *hTrkVertexY = 
455     new TH1F("hTrkVertexY","ITS+TPC Trk Vertex y; y [cm]; N events",
456              10000,-2,2);
457   hTrkVertexY->Sumw2();
458   Add2ESDsList(hTrkVertexY, 9, expertHistogram);
459
460   TH1F *hTrkVertexZ = 
461     new TH1F("hTrkVertexZ","ITS+TPC Trk Vertex Z; z [cm]; N events",
462              10000,-20,20);
463   hTrkVertexZ->Sumw2();
464   Add2ESDsList(hTrkVertexZ, 10, expertHistogram);
465
466   TH1F *hTrkVertexContrOverITSrefit5 =
467     new TH1F("hTrkVertexContrOverITSrefit5","ITS+TPC Trk Vertex: contributors / tracks; N contributors / N trks kITSrefit with 5 or 6 clusters; N events",
468              100,-4,2);
469   hTrkVertexContrOverITSrefit5->Sumw2();
470   Add2ESDsList(hTrkVertexContrOverITSrefit5, 11, expertHistogram);
471
472   TH1F *hSPDTrkVertexDeltaX =
473     new TH1F("hSPDTrkVertexDeltaX","Comparison of SPD and Trk vertices: x; xSPD-xTrk [cm]; N events",
474              1000,-1,1);
475   hSPDTrkVertexDeltaX->Sumw2();
476   Add2ESDsList(hSPDTrkVertexDeltaX, 12, expertHistogram);
477     
478   TH1F *hSPDTrkVertexDeltaY =
479     new TH1F("hSPDTrkVertexDeltaY","Comparison of SPD and Trk vertices: y; ySPD-yTrk [cm]; N events",
480              1000,-1,1);
481   hSPDTrkVertexDeltaY->Sumw2();
482   Add2ESDsList(hSPDTrkVertexDeltaY, 13, expertHistogram);
483     
484   TH1F *hSPDTrkVertexDeltaZ =
485     new TH1F("hSPDTrkVertexDeltaZ","Comparison of SPD and Trk vertices: z; zSPD-zTrk [cm]; N events",
486              1000,-1,1);
487   hSPDTrkVertexDeltaZ->Sumw2();
488   Add2ESDsList(hSPDTrkVertexDeltaZ, 14);
489     
490   // SPD Tracklets
491
492   TH1F* hSPDTracklets = 
493     new TH1F("hSPDTracklets","N SPD Tracklets; N tracklets; Counts",300,0.,300.);
494   hSPDTracklets->Sumw2();
495   Add2ESDsList(hSPDTracklets, 15); 
496
497   TH2F* hSPDTrackletsvsFiredChips0 = 
498     new TH2F("hSPDTrackletsvsFiredChips0","N SPD Tracklets vs N FiredChips Layer0",
499               300,0.,300.,300,0.,300.);
500   hSPDTrackletsvsFiredChips0->GetXaxis()->SetTitle("N FiredChips Layer0"); 
501   hSPDTrackletsvsFiredChips0->GetYaxis()->SetTitle("N SPD Tracklets"); 
502   hSPDTrackletsvsFiredChips0->Sumw2();
503   Add2ESDsList(hSPDTrackletsvsFiredChips0, 16, expertHistogram ); 
504
505   TH2F* hSPDTrackletsvsFiredChips1 = 
506     new TH2F("hSPDTrackletsvsFiredChips1","N SPD Tracklets vs N FiredChips Layer1",
507               300,0.,300.,300,0.,300.);
508   hSPDTrackletsvsFiredChips1->GetXaxis()->SetTitle("N FiredChips Layer1"); 
509   hSPDTrackletsvsFiredChips1->GetYaxis()->SetTitle("N SPD Tracklets"); 
510   hSPDTrackletsvsFiredChips1->Sumw2();
511   Add2ESDsList(hSPDTrackletsvsFiredChips1, 17, expertHistogram); 
512
513   TH2F* hSPDFiredChips1vsFiredChips0 = 
514     new TH2F("hSPDFiredChips1vsFiredChips0","N FiredChips Layer1 vs N FiredChips Layer0",
515               300,0.,300.,300,0.,300.);
516   hSPDFiredChips1vsFiredChips0->GetXaxis()->SetTitle("N FiredChips Layer0"); 
517   hSPDFiredChips1vsFiredChips0->GetYaxis()->SetTitle("N FiredChips Layer1"); 
518   hSPDFiredChips1vsFiredChips0->Sumw2();
519   Add2ESDsList(hSPDFiredChips1vsFiredChips0, 18, expertHistogram ); 
520     
521   TH1F* hSPDTrackletsDePhi = 
522     new TH1F("hSPDTrackletsDePhi","DeltaPhi SPD Tracklets; DeltaPhi [rad]; N events",200,-0.2,0.2);
523   hSPDTrackletsDePhi->Sumw2();
524   Add2ESDsList(hSPDTrackletsDePhi, 19); 
525     
526   TH1F* hSPDTrackletsPhi = 
527     new TH1F("hSPDTrackletsPhi","Phi SPD Tracklets; Phi [rad]; N events",1000,0.,2*TMath::Pi());
528   hSPDTrackletsPhi->Sumw2();
529   Add2ESDsList(hSPDTrackletsPhi, 20); 
530     
531   TH1F* hSPDTrackletsDeTheta = 
532     new TH1F("hSPDTrackletsDeTheta","DeltaTheta SPD Tracklets; DeltaTheta [rad]; N events",200,-0.2,0.2);
533   hSPDTrackletsDeTheta->Sumw2();
534   Add2ESDsList(hSPDTrackletsDeTheta, 21); 
535
536   TH1F* hSPDTrackletsTheta = 
537     new TH1F("hSPDTrackletsTheta","Theta SPD Tracklets; Theta [rad]; N events",500,0.,TMath::Pi());
538   hSPDTrackletsTheta->Sumw2();
539   Add2ESDsList(hSPDTrackletsTheta, 22); 
540
541   // map of layers skipped by tracking (set in AliITSRecoParam)
542   TH1F *hESDSkippedLayers = 
543     new TH1F("hESDSkippedLayers", "Map of layers skipped by tracking; Layer; Skipped",
544              6, -0.5, 5.5);
545   hESDSkippedLayers->Sumw2();
546   hESDSkippedLayers->SetMinimum(0);
547   Add2ESDsList(hESDSkippedLayers, 23, expertHistogram);
548
549
550   return;
551 }
552
553 //____________________________________________________________________________
554 void AliITSQADataMakerRec::MakeESDs(AliESDEvent *esd)
555 {
556   // Make QA data from ESDs
557
558   // Check id histograms already created for this Event Specie
559 //  if ( ! GetESDsData(0) )
560 //    InitESDs() ;
561   
562   const Int_t nESDTracks = esd->GetNumberOfTracks();
563   Int_t nITSrefit5 = 0; 
564
565   Int_t idet,status;
566   Float_t xloc,zloc;
567
568   // loop on tracks
569   for(Int_t i = 0; i < nESDTracks; i++) {
570     
571     AliESDtrack *track = esd->GetTrack(i);
572     
573     Int_t nclsITS = track->GetNcls(0);
574
575     Bool_t itsrefit=kFALSE,tpcin=kFALSE,itsin=kFALSE;
576     if ((track->GetStatus() & AliESDtrack::kITSrefit)) itsrefit=kTRUE;
577     if ((track->GetStatus() & AliESDtrack::kTPCin)) tpcin=kTRUE;     
578     if ((track->GetStatus() & AliESDtrack::kITSin)) itsin=kTRUE;     
579     if(nclsITS>=5 && itsrefit) nITSrefit5++;
580
581     if(tpcin) {
582       GetESDsData(0)->Fill(nclsITS);
583     }
584     if(itsin && !tpcin){
585       GetESDsData(2)->Fill(nclsITS);
586     }
587
588     for(Int_t layer=0; layer<6; layer++) {
589
590       if(TESTBIT(track->GetITSClusterMap(),layer)) {
591         if(tpcin) {
592           GetESDsData(1)->Fill(layer);
593         } else {
594           GetESDsData(3)->Fill(layer);
595         }
596       }
597       track->GetITSModuleIndexInfo(layer,idet,status,xloc,zloc);
598       if(status==3) GetESDsData(23)->SetBinContent(layer,1);
599     }     
600
601   } // end loop on tracks
602
603   // vertices
604   const AliESDVertex *vtxSPD = esd->GetPrimaryVertexSPD();
605   const AliESDVertex *vtxTrk = esd->GetPrimaryVertexTracks();
606
607   Int_t mult = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetNumberOfTracklets();
608   if(mult>0)
609     GetESDsData(7)->Fill((Float_t)(vtxSPD->GetNContributors())/(Float_t)mult);
610
611   if(nITSrefit5>0)
612     GetESDsData(11)->Fill((Float_t)(vtxTrk->GetNIndices())/(Float_t)nITSrefit5);
613
614   if(vtxSPD->GetNContributors()>0) {
615     GetESDsData(4)->Fill(vtxSPD->GetXv());
616     GetESDsData(5)->Fill(vtxSPD->GetYv());
617     GetESDsData(6)->Fill(vtxSPD->GetZv());
618   }
619
620   if(vtxTrk->GetNContributors()>0) {
621     GetESDsData(8)->Fill(vtxTrk->GetXv());
622     GetESDsData(9)->Fill(vtxTrk->GetYv());
623     GetESDsData(10)->Fill(vtxTrk->GetZv());
624   }
625
626   if(vtxSPD->GetNContributors()>0 && 
627      vtxTrk->GetNContributors()>0) {
628     GetESDsData(12)->Fill(vtxSPD->GetXv()-vtxTrk->GetXv());
629     GetESDsData(13)->Fill(vtxSPD->GetYv()-vtxTrk->GetYv());
630     GetESDsData(14)->Fill(vtxSPD->GetZv()-vtxTrk->GetZv());
631   }
632
633   // SPD Tracklets
634   GetESDsData(15)->Fill(mult);
635
636   Short_t nFiredChips0 = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetNumberOfFiredChips(0);
637   Short_t nFiredChips1 = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetNumberOfFiredChips(1);
638   GetESDsData(16)->Fill(nFiredChips0,mult);
639   GetESDsData(17)->Fill(nFiredChips1,mult);
640   GetESDsData(18)->Fill(nFiredChips0,nFiredChips1);
641
642   // Loop over tracklets
643   for (Int_t itr=0; itr<mult; ++itr) {
644     Float_t dePhiTr   = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetDeltaPhi(itr);
645     Float_t deThetaTr = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetDeltaTheta(itr);
646     Float_t phiTr   = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetPhi(itr);
647     Float_t thetaTr = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetTheta(itr);
648     GetESDsData(19)->Fill(dePhiTr);
649     GetESDsData(20)->Fill(phiTr);
650     GetESDsData(21)->Fill(deThetaTr);
651     GetESDsData(22)->Fill(thetaTr);
652   } // end loop on tracklets
653
654   return;
655 }
656
657 //_________________________________________________________________
658 Int_t AliITSQADataMakerRec::GetDetTaskOffset(Int_t subdet,AliQAv1::TASKINDEX_t task)
659 {
660   switch(subdet)
661     {
662
663       Int_t offset;
664     case 1:
665       offset=fSPDDataMaker->GetOffset(task);
666       return offset;
667       break;
668     case 2:
669       offset=fSDDDataMaker->GetOffset(task);
670       return offset;
671       break;
672     case 3:
673       offset=fSSDDataMaker->GetOffset(task);
674       return offset;
675       break;
676     default:
677       AliWarning("No specific subdetector (SPD, SDD, SSD) selected!! Offset set to zero \n");
678       offset=0;
679       return offset;
680       break;
681     }
682   //return offset;
683 }
684
685 //____________________________________________________________________
686
687 Bool_t AliITSQADataMakerRec::AreEqual(Double_t a1,Double_t a2)
688 {
689   const Double_t kEpsilon= 1.e-14;
690   return TMath::Abs(a1-a2)<=kEpsilon*TMath::Abs(a1);      
691 }
692