]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSQADataMakerRec.cxx
Updated version of ITS QA Checker and related modifications (Melinda)
[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   AliInfo(Form("End of Dedetctor Cycle called for %s\n",AliQAv1::GetTaskName(task).Data() ));
134   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
135
136         if(AliQAv1::Instance()->IsEventSpecieSet(specie)){
137           Int_t idnumber=list[specie]->GetUniqueID();
138           //printf("specie %s \t id number == %d\n",AliRecoParam::GetEventSpecieName(specie),idnumber);
139           if(idnumber==40||idnumber==0)
140             {
141               //AliInfo(Form("No check for %s\n",AliQAv1::GetTaskName(task).Data() ))
142                 continue;
143             } //skip kDigitsR and not filled TobjArray specie
144           else{
145             AliDebug(AliQAv1::GetQADebugLevel(),"AliITSDM instantiates checker with Run(AliQAv1::kITS, task, list[specie])\n"); 
146             if(fSubDetector == 0 || fSubDetector == 1) fSPDDataMaker->EndOfDetectorCycle(task, list[/*GetEventSpecie()*/specie]);
147             if(fSubDetector == 0 || fSubDetector == 2) fSDDDataMaker->EndOfDetectorCycle(task, list[/*GetEventSpecie()*/specie]);
148             if(fSubDetector == 0 || fSubDetector == 3) fSSDDataMaker->EndOfDetectorCycle(task, list[/*GetEventSpecie()*/specie]);
149             
150             
151             AliQAChecker *qac = AliQAChecker::Instance();
152             AliITSQAChecker *qacb = (AliITSQAChecker *) qac->GetDetQAChecker(0);
153             Int_t subdet=GetSubDet();
154             qacb->SetSubDet(subdet);
155           
156             if(subdet== 0 ){
157               qacb->SetTaskOffset(fSPDDataMaker->GetOffset(task,specie), fSDDDataMaker->GetOffset(task,specie), fSSDDataMaker->GetOffset(task,specie)); //Setting the offset for the QAChecker list
158               qacb->SetHisto(fSPDDataMaker->GetTaskHisto(task), fSDDDataMaker->GetTaskHisto(task), fSSDDataMaker->GetTaskHisto(task));
159             }
160             else
161               if(subdet!=0){
162                 Int_t offset=GetDetTaskOffset(subdet, task,specie);
163                 qacb->SetDetTaskOffset(subdet,offset);
164                 Int_t histo=GetDetTaskHisto(subdet, task);
165                 qacb->SetDetHisto(subdet,histo);
166               }
167             
168             qac->Run( AliQAv1::kITS , task, list);
169             
170           }//end else unique id
171         }//end else event specie
172   }//end for
173 }
174
175 //____________________________________________________________________________ 
176 void AliITSQADataMakerRec::EndOfDetectorCycle(const char * /*fgDataName*/)
177 {
178   //eventually used for different  AliQAChecker::Instance()->Run
179 }
180
181 //____________________________________________________________________________ 
182 void AliITSQADataMakerRec::InitRaws()
183 {  
184
185   //if(fRawsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries()) return;
186         
187
188
189   if(fSubDetector == 0 || fSubDetector == 1) {
190     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SPD InitRaws\n");
191     fSPDDataMaker->InitRaws();
192   }
193   if(fSubDetector == 0 || fSubDetector == 2) {
194     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SDD InitRaws\n");
195     
196     fSDDDataMaker->SetOffset(AliQAv1::kRAWS, fRawsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries(),AliRecoParam::AConvert(fEventSpecie));
197     fSDDDataMaker->InitRaws();
198   }
199   if(fSubDetector == 0 || fSubDetector == 3) {
200     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SSD InitRaws\n");
201     
202     fSSDDataMaker->SetOffset(AliQAv1::kRAWS, fRawsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries(),AliRecoParam::AConvert(fEventSpecie));
203     fSSDDataMaker->InitRaws();
204   }
205   fRawsQAList[AliRecoParam::AConvert(fEventSpecie)]->SetUniqueID(10);
206
207 }
208
209 //____________________________________________________________________________
210 void AliITSQADataMakerRec::MakeRaws(AliRawReader* rawReader)
211
212   // Fill QA for RAW   
213   //return ; 
214
215   SetRunNumber(rawReader->GetRunNumber());
216
217   if(fSubDetector == 0 || fSubDetector == 1)  {
218     fSPDDataMaker->MakeRaws(rawReader) ; 
219   }
220   
221   if(fSubDetector == 0 || fSubDetector == 2) {
222     fSDDDataMaker->MakeRaws(rawReader) ; 
223   }
224
225   if(fSubDetector == 0 || fSubDetector == 3) fSSDDataMaker->MakeRaws(rawReader);
226
227 }
228
229 //____________________________________________________________________________ 
230 void AliITSQADataMakerRec::InitDigits()
231 {
232
233   // Initialization for DIGITS
234   if(fSubDetector == 0 || fSubDetector == 1) {
235     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SPD InitDigits\n");
236
237     fSPDDataMaker->InitDigits();
238   }
239   if(fSubDetector == 0 || fSubDetector == 2) {
240     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SDD InitDigits\n");
241     fSDDDataMaker->SetOffset(AliQAv1::kDIGITSR, fDigitsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries(),AliRecoParam::AConvert(fEventSpecie));
242
243     fSDDDataMaker->InitDigits();
244   }
245   if(fSubDetector == 0 || fSubDetector == 3) {
246     AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SSD InitDigits\n");
247     fSSDDataMaker->SetOffset(AliQAv1::kDIGITSR, fDigitsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries(),AliRecoParam::AConvert(fEventSpecie));
248
249     fSSDDataMaker->InitDigits();
250   }
251   fDigitsQAList[AliRecoParam::AConvert(fEventSpecie)]->SetUniqueID(40);
252 }
253
254 //____________________________________________________________________________ 
255 void AliITSQADataMakerRec::MakeDigits(TTree * digitsTree)
256 {
257
258   
259   // Fill QA for recpoints
260   if(fSubDetector == 0 || fSubDetector == 1) {
261     fSPDDataMaker->MakeDigits(digitsTree) ; 
262   }
263   
264   if(fSubDetector == 0 || fSubDetector == 2) {
265     fSDDDataMaker->MakeDigits(digitsTree) ; 
266     
267   }
268   
269   if(fSubDetector == 0 || fSubDetector == 3) fSSDDataMaker->MakeDigits(digitsTree);
270 }
271
272 //____________________________________________________________________________ 
273 void AliITSQADataMakerRec::InitRecPoints()
274 {
275
276   // Initialization for RECPOINTS
277
278
279   //if(fRecPointsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries()) return;
280         if(fSubDetector == 0 || fSubDetector == 1) {
281                 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SPD InitRecPoints\n");
282                 fSPDDataMaker->InitRecPoints();
283         }
284         if(fSubDetector == 0 || fSubDetector == 2) {
285                 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SDD InitRecPoints\n");
286                 fSDDDataMaker->SetOffset(AliQAv1::kRECPOINTS, fRecPointsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries(), AliRecoParam::AConvert(fEventSpecie));
287                 fSDDDataMaker->InitRecPoints();
288         }
289         if(fSubDetector == 0 || fSubDetector == 3) {
290                 AliDebug(AliQAv1::GetQADebugLevel(),"AliITSQADM:: SSD InitRecPoints\n");
291                 fSSDDataMaker->SetOffset(AliQAv1::kRECPOINTS, fRecPointsQAList[AliRecoParam::AConvert(fEventSpecie)]->GetEntries(),AliRecoParam::AConvert(fEventSpecie));
292                 fSSDDataMaker->InitRecPoints();
293         }
294
295   fRecPointsQAList[AliRecoParam::AConvert(fEventSpecie)]->SetUniqueID(20);
296         if(fSubDetector == 0){
297           Int_t offset = fRecPointsQAList [AliRecoParam::AConvert(fEventSpecie)]->GetEntries();
298           const Bool_t expert   = kTRUE ; 
299           const Bool_t image    = kTRUE ; 
300           Char_t name[50];
301           Char_t title[50];
302           TH2F**hPhiEta = new TH2F*[6];
303           for (Int_t iLay=0;iLay<6;iLay++) {
304             sprintf(name,"Phi_vs_Eta_ITS_Layer%d",iLay+1);
305             sprintf(title,"Phi vs Eta - ITS Layer %d",iLay+1);
306             hPhiEta[iLay]=new TH2F(name,title,30,-1.5,1.5,200,0.,2*TMath::Pi());
307             hPhiEta[iLay]->GetXaxis()->SetTitle("Pseudorapidity");
308             hPhiEta[iLay]->GetYaxis()->SetTitle("#varphi [rad]");
309             Add2RecPointsList(hPhiEta[iLay], iLay + offset, !expert, image);
310             
311             //delete hPhiEta[iLay];
312           }
313           
314         }
315         
316 }
317
318 //____________________________________________________________________________ 
319 void AliITSQADataMakerRec::MakeRecPoints(TTree * clustersTree)
320 {
321  
322   // Fill QA for recpoints
323   if(fSubDetector == 0 || fSubDetector == 1) {
324     fSPDDataMaker->MakeRecPoints(clustersTree) ; 
325   }
326     
327   if(fSubDetector == 0 || fSubDetector == 2) {
328     fSDDDataMaker->MakeRecPoints(clustersTree) ; 
329   }
330   
331   if(fSubDetector == 0 || fSubDetector == 3) fSSDDataMaker->MakeRecPoints(clustersTree);
332
333
334
335   if(fSubDetector == 0){
336     // Check id histograms already created for this Event Specie
337     TBranch *branchRecP = clustersTree->GetBranch("ITSRecPoints");
338     if (!branchRecP) {
339       AliError("can't get the branch with the ITS clusters !");
340       return;
341     }
342     
343     Int_t offset = fRecPointsQAList [AliRecoParam::AConvert(fEventSpecie)]->GetEntries();
344     Float_t cluGlo[3] = {0.,0.,0.};
345     Int_t lay, lad, det; 
346     static TClonesArray statRecpoints("AliITSRecPoint") ;
347     TClonesArray *recpoints = &statRecpoints;
348     branchRecP->SetAddress(&recpoints);
349     // Fill QA for recpoints
350     for(Int_t module=0; module<clustersTree->GetEntries();module++){
351       //  AliInfo(Form("Module %d\n",module));
352       branchRecP->GetEvent(module);
353       AliITSgeomTGeo::GetModuleId(module, lay, lad, det);
354       for(Int_t j=0;j<recpoints->GetEntries();j++){
355         AliITSRecPoint *rcp = (AliITSRecPoint*)recpoints->At(j);    
356         //Check id histograms already created for this Event Specie
357         rcp->GetGlobalXYZ(cluGlo);
358         Double_t rad=TMath::Sqrt(cluGlo[0]*cluGlo[0]+cluGlo[1]*cluGlo[1]+cluGlo[2]*cluGlo[2]);
359         Double_t phi= TMath::Pi() + TMath::ATan2(-cluGlo[1],-cluGlo[0]);
360         Double_t theta = TMath::ACos(cluGlo[2]/rad);
361         Double_t eta = 100.;
362         if(AreEqual(rad,0.) == kFALSE) {
363           if(theta<=1.e-14){ eta=30.; }
364           else { eta = -TMath::Log(TMath::Tan(theta/2.));}
365         }
366         //      printf("=========================>hlt   rcp->GetLayer() = %d \n",rcp->GetLayer());
367         (GetRecPointsData( rcp->GetLayer() + offset - 6))->Fill(eta,phi);
368       }
369     }
370   }
371   
372 }
373
374 //____________________________________________________________________________ 
375 void AliITSQADataMakerRec::FillRecPoint(AliITSRecPoint rcp)
376 {
377
378         // Fill QA for recpoints
379         Float_t cluGlo[3] = {0.,0.,0.};
380         Int_t offset = fRecPointsQAList [AliRecoParam::AConvert(fEventSpecie)]->GetEntries();
381   // Check id histograms already created for this Event Specie
382         rcp.GetGlobalXYZ(cluGlo);
383         Double_t rad=TMath::Sqrt(cluGlo[0]*cluGlo[0]+cluGlo[1]*cluGlo[1]+cluGlo[2]*cluGlo[2]);
384         Double_t phi= TMath::Pi() + TMath::ATan2(-cluGlo[1],-cluGlo[0]);
385         Double_t theta = TMath::ACos(cluGlo[2]/rad);
386         Double_t eta = 100.;
387         if(AreEqual(rad,0.)==kFALSE) {
388           if(theta<=1.e-14){eta=30.;}
389           else    {eta = -TMath::Log(TMath::Tan(theta/2.));}
390         }
391         (GetRecPointsData( rcp.GetLayer() + offset - 6))->Fill(eta,phi);        
392
393 }
394
395 //____________________________________________________________________________ 
396 TH2F *AliITSQADataMakerRec::GetITSGlobalHisto(Int_t layer)
397 {
398
399         Int_t offset = fRecPointsQAList [AliRecoParam::AConvert(fEventSpecie)]->GetEntries();
400         return ((TH2F *) GetRecPointsData( layer + offset - 6));//local distribution
401 }
402
403 //____________________________________________________________________________ 
404 void AliITSQADataMakerRec::InitESDs()
405 {
406
407   // Create ESDs histograms in ESDs subdir
408
409   Bool_t expertHistogram = kTRUE;
410
411   TH1F *hESDClustersMI = 
412     new TH1F("hESDClustersMI", "N ITS clusters per track (MI); N clusters; Counts",
413              7, -0.5, 6.5);
414   hESDClustersMI->Sumw2();
415   hESDClustersMI->SetMinimum(0);
416   Add2ESDsList(hESDClustersMI, 0);
417
418   TH1F *hESDClusterMapMI =
419     new TH1F("hESDClusterMapMI", "N tracks with point Layer (MI); Layer; N tracks",
420              6, -0.5, 5.5);
421   hESDClusterMapMI->Sumw2();
422   hESDClusterMapMI->SetMinimum(0);
423   Add2ESDsList(hESDClusterMapMI, 1, expertHistogram);
424
425   TH1F *hESDClustersSA = 
426     new TH1F("hESDClustersSA", "N ITS clusters per track (SA); N clusters; Counts",
427              7, -0.5, 6.5);
428   hESDClustersSA->Sumw2();
429   hESDClustersSA->SetMinimum(0);
430   Add2ESDsList(hESDClustersSA, 2);
431
432   TH1F *hESDClusterMapSA =
433     new TH1F("hESDClusterMapSA", "N tracks with point Layer (SA); Layer; N tracks",
434              6, -0.5, 5.5);
435   hESDClusterMapSA->Sumw2();
436   hESDClusterMapSA->SetMinimum(0);
437   Add2ESDsList(hESDClusterMapSA, 3, expertHistogram);
438
439   TH1F *hSPDVertexX = 
440     new TH1F("hSPDVertexX","SPD Vertex x; x [cm]; N events",
441              10000,-2,2);
442   hSPDVertexX->Sumw2();
443   Add2ESDsList(hSPDVertexX, 4);
444
445   TH1F *hSPDVertexY = 
446     new TH1F("hSPDVertexY","SPD Vertex y; y [cm]; N events",
447              10000,-2,2);
448   hSPDVertexY->Sumw2();
449   Add2ESDsList(hSPDVertexY, 5);
450
451   TH1F *hSPDVertexZ = 
452     new TH1F("hSPDVertexZ","SPD Vertex Z; z [cm]; N events",
453              10000,-20,20);
454   hSPDVertexZ->Sumw2();
455   Add2ESDsList(hSPDVertexZ, 6);
456
457   TH1F *hSPDVertexContrOverMult =
458     new TH1F("hSPDVertexContrOverMult","SPD Vertex: contributors / multiplicity; N contributors / SPD multiplicity; N events",
459              100,-4,20);
460   hSPDVertexContrOverMult->Sumw2();
461   Add2ESDsList(hSPDVertexContrOverMult, 7, expertHistogram);
462
463   TH1F *hTrkVertexX = 
464     new TH1F("hTrkVertexX","ITS+TPC Trk Vertex x; x [cm]; N events",
465              10000,-2,2);
466   hTrkVertexX->Sumw2();
467   Add2ESDsList(hTrkVertexX, 8, expertHistogram);
468
469   TH1F *hTrkVertexY = 
470     new TH1F("hTrkVertexY","ITS+TPC Trk Vertex y; y [cm]; N events",
471              10000,-2,2);
472   hTrkVertexY->Sumw2();
473   Add2ESDsList(hTrkVertexY, 9, expertHistogram);
474
475   TH1F *hTrkVertexZ = 
476     new TH1F("hTrkVertexZ","ITS+TPC Trk Vertex Z; z [cm]; N events",
477              10000,-20,20);
478   hTrkVertexZ->Sumw2();
479   Add2ESDsList(hTrkVertexZ, 10, expertHistogram);
480
481   TH1F *hTrkVertexContrOverITSrefit5 =
482     new TH1F("hTrkVertexContrOverITSrefit5","ITS+TPC Trk Vertex: contributors / tracks; N contributors / N trks kITSrefit with 5 or 6 clusters; N events",
483              100,-4,2);
484   hTrkVertexContrOverITSrefit5->Sumw2();
485   Add2ESDsList(hTrkVertexContrOverITSrefit5, 11, expertHistogram);
486
487   TH1F *hSPDTrkVertexDeltaX =
488     new TH1F("hSPDTrkVertexDeltaX","Comparison of SPD and Trk vertices: x; xSPD-xTrk [cm]; N events",
489              1000,-1,1);
490   hSPDTrkVertexDeltaX->Sumw2();
491   Add2ESDsList(hSPDTrkVertexDeltaX, 12, expertHistogram);
492     
493   TH1F *hSPDTrkVertexDeltaY =
494     new TH1F("hSPDTrkVertexDeltaY","Comparison of SPD and Trk vertices: y; ySPD-yTrk [cm]; N events",
495              1000,-1,1);
496   hSPDTrkVertexDeltaY->Sumw2();
497   Add2ESDsList(hSPDTrkVertexDeltaY, 13, expertHistogram);
498     
499   TH1F *hSPDTrkVertexDeltaZ =
500     new TH1F("hSPDTrkVertexDeltaZ","Comparison of SPD and Trk vertices: z; zSPD-zTrk [cm]; N events",
501              1000,-1,1);
502   hSPDTrkVertexDeltaZ->Sumw2();
503   Add2ESDsList(hSPDTrkVertexDeltaZ, 14);
504     
505   // SPD Tracklets
506
507   TH1F* hSPDTracklets = 
508     new TH1F("hSPDTracklets","N SPD Tracklets; N tracklets; Counts",300,0.,300.);
509   hSPDTracklets->Sumw2();
510   Add2ESDsList(hSPDTracklets, 15); 
511
512   TH2F* hSPDTrackletsvsFiredChips0 = 
513     new TH2F("hSPDTrackletsvsFiredChips0","N SPD Tracklets vs N FiredChips Layer0",
514               300,0.,300.,300,0.,300.);
515   hSPDTrackletsvsFiredChips0->GetXaxis()->SetTitle("N FiredChips Layer0"); 
516   hSPDTrackletsvsFiredChips0->GetYaxis()->SetTitle("N SPD Tracklets"); 
517   hSPDTrackletsvsFiredChips0->Sumw2();
518   Add2ESDsList(hSPDTrackletsvsFiredChips0, 16, expertHistogram ); 
519
520   TH2F* hSPDTrackletsvsFiredChips1 = 
521     new TH2F("hSPDTrackletsvsFiredChips1","N SPD Tracklets vs N FiredChips Layer1",
522               300,0.,300.,300,0.,300.);
523   hSPDTrackletsvsFiredChips1->GetXaxis()->SetTitle("N FiredChips Layer1"); 
524   hSPDTrackletsvsFiredChips1->GetYaxis()->SetTitle("N SPD Tracklets"); 
525   hSPDTrackletsvsFiredChips1->Sumw2();
526   Add2ESDsList(hSPDTrackletsvsFiredChips1, 17, expertHistogram); 
527
528   TH2F* hSPDFiredChips1vsFiredChips0 = 
529     new TH2F("hSPDFiredChips1vsFiredChips0","N FiredChips Layer1 vs N FiredChips Layer0",
530               300,0.,300.,300,0.,300.);
531   hSPDFiredChips1vsFiredChips0->GetXaxis()->SetTitle("N FiredChips Layer0"); 
532   hSPDFiredChips1vsFiredChips0->GetYaxis()->SetTitle("N FiredChips Layer1"); 
533   hSPDFiredChips1vsFiredChips0->Sumw2();
534   Add2ESDsList(hSPDFiredChips1vsFiredChips0, 18, expertHistogram ); 
535     
536   TH1F* hSPDTrackletsDePhi = 
537     new TH1F("hSPDTrackletsDePhi","DeltaPhi SPD Tracklets; DeltaPhi [rad]; N events",200,-0.2,0.2);
538   hSPDTrackletsDePhi->Sumw2();
539   Add2ESDsList(hSPDTrackletsDePhi, 19); 
540     
541   TH1F* hSPDTrackletsPhi = 
542     new TH1F("hSPDTrackletsPhi","Phi SPD Tracklets; Phi [rad]; N events",1000,0.,2*TMath::Pi());
543   hSPDTrackletsPhi->Sumw2();
544   Add2ESDsList(hSPDTrackletsPhi, 20); 
545     
546   TH1F* hSPDTrackletsDeTheta = 
547     new TH1F("hSPDTrackletsDeTheta","DeltaTheta SPD Tracklets; DeltaTheta [rad]; N events",200,-0.2,0.2);
548   hSPDTrackletsDeTheta->Sumw2();
549   Add2ESDsList(hSPDTrackletsDeTheta, 21); 
550
551   TH1F* hSPDTrackletsTheta = 
552     new TH1F("hSPDTrackletsTheta","Theta SPD Tracklets; Theta [rad]; N events",500,0.,TMath::Pi());
553   hSPDTrackletsTheta->Sumw2();
554   Add2ESDsList(hSPDTrackletsTheta, 22); 
555
556   // map of layers skipped by tracking (set in AliITSRecoParam)
557   TH1F *hESDSkippedLayers = 
558     new TH1F("hESDSkippedLayers", "Map of layers skipped by tracking; Layer; Skipped",
559              6, -0.5, 5.5);
560   hESDSkippedLayers->Sumw2();
561   hESDSkippedLayers->SetMinimum(0);
562   Add2ESDsList(hESDSkippedLayers, 23, expertHistogram);
563
564   fESDsQAList[AliRecoParam::AConvert(fEventSpecie)]->SetUniqueID(30);
565   return;
566 }
567
568 //____________________________________________________________________________
569 void AliITSQADataMakerRec::MakeESDs(AliESDEvent *esd)
570 {
571   // Make QA data from ESDs
572
573   // Check id histograms already created for this Event Specie
574 //  if ( ! GetESDsData(0) )
575 //    InitESDs() ;
576   
577   const Int_t nESDTracks = esd->GetNumberOfTracks();
578   Int_t nITSrefit5 = 0; 
579
580   Int_t idet,status;
581   Float_t xloc,zloc;
582
583   // loop on tracks
584   for(Int_t i = 0; i < nESDTracks; i++) {
585     
586     AliESDtrack *track = esd->GetTrack(i);
587     
588     Int_t nclsITS = track->GetNcls(0);
589
590     Bool_t itsrefit=kFALSE,tpcin=kFALSE,itsin=kFALSE;
591     if ((track->GetStatus() & AliESDtrack::kITSrefit)) itsrefit=kTRUE;
592     if ((track->GetStatus() & AliESDtrack::kTPCin)) tpcin=kTRUE;     
593     if ((track->GetStatus() & AliESDtrack::kITSin)) itsin=kTRUE;     
594     if(nclsITS>=5 && itsrefit) nITSrefit5++;
595
596     if(tpcin) {
597       GetESDsData(0)->Fill(nclsITS);
598     }
599     if(itsin && !tpcin){
600       GetESDsData(2)->Fill(nclsITS);
601     }
602
603     for(Int_t layer=0; layer<6; layer++) {
604
605       if(TESTBIT(track->GetITSClusterMap(),layer)) {
606         if(tpcin) {
607           GetESDsData(1)->Fill(layer);
608         } else {
609           GetESDsData(3)->Fill(layer);
610         }
611       }
612       track->GetITSModuleIndexInfo(layer,idet,status,xloc,zloc);
613       if(status==3) GetESDsData(23)->SetBinContent(layer,1);
614     }     
615
616   } // end loop on tracks
617
618   // vertices
619   const AliESDVertex *vtxSPD = esd->GetPrimaryVertexSPD();
620   const AliESDVertex *vtxTrk = esd->GetPrimaryVertexTracks();
621
622   Int_t mult = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetNumberOfTracklets();
623   if(mult>0)
624     GetESDsData(7)->Fill((Float_t)(vtxSPD->GetNContributors())/(Float_t)mult);
625
626   if(nITSrefit5>0)
627     GetESDsData(11)->Fill((Float_t)(vtxTrk->GetNIndices())/(Float_t)nITSrefit5);
628
629   if(vtxSPD->GetNContributors()>0) {
630     GetESDsData(4)->Fill(vtxSPD->GetXv());
631     GetESDsData(5)->Fill(vtxSPD->GetYv());
632     GetESDsData(6)->Fill(vtxSPD->GetZv());
633   }
634
635   if(vtxTrk->GetNContributors()>0) {
636     GetESDsData(8)->Fill(vtxTrk->GetXv());
637     GetESDsData(9)->Fill(vtxTrk->GetYv());
638     GetESDsData(10)->Fill(vtxTrk->GetZv());
639   }
640
641   if(vtxSPD->GetNContributors()>0 && 
642      vtxTrk->GetNContributors()>0) {
643     GetESDsData(12)->Fill(vtxSPD->GetXv()-vtxTrk->GetXv());
644     GetESDsData(13)->Fill(vtxSPD->GetYv()-vtxTrk->GetYv());
645     GetESDsData(14)->Fill(vtxSPD->GetZv()-vtxTrk->GetZv());
646   }
647
648   // SPD Tracklets
649   GetESDsData(15)->Fill(mult);
650
651   Short_t nFiredChips0 = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetNumberOfFiredChips(0);
652   Short_t nFiredChips1 = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetNumberOfFiredChips(1);
653   GetESDsData(16)->Fill(nFiredChips0,mult);
654   GetESDsData(17)->Fill(nFiredChips1,mult);
655   GetESDsData(18)->Fill(nFiredChips0,nFiredChips1);
656
657   // Loop over tracklets
658   for (Int_t itr=0; itr<mult; ++itr) {
659     Float_t dePhiTr   = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetDeltaPhi(itr);
660     Float_t deThetaTr = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetDeltaTheta(itr);
661     Float_t phiTr   = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetPhi(itr);
662     Float_t thetaTr = ((AliMultiplicity*)(esd->GetMultiplicity()))->GetTheta(itr);
663     GetESDsData(19)->Fill(dePhiTr);
664     GetESDsData(20)->Fill(phiTr);
665     GetESDsData(21)->Fill(deThetaTr);
666     GetESDsData(22)->Fill(thetaTr);
667   } // end loop on tracklets
668
669   return;
670 }
671
672 //_________________________________________________________________
673 Int_t AliITSQADataMakerRec::GetDetTaskOffset(Int_t subdet,AliQAv1::TASKINDEX_t task, Int_t specie)
674 {
675
676   Int_t offset=0;
677   switch(subdet)
678     {
679     case 1:
680       offset=fSPDDataMaker->GetOffset(task,specie);
681       //return offset;
682       break;
683     case 2:
684       offset=fSDDDataMaker->GetOffset(task,specie);
685       //return offset;
686       break;
687     case 3:
688       offset=fSSDDataMaker->GetOffset(task,specie);
689       //return offset;
690       break;
691     default:
692       AliWarning("No specific subdetector (SPD, SDD, SSD) selected!! Offset set to zero \n");
693       offset=0;
694       //return offset;
695       break;
696     }
697   return offset;
698 }
699
700 //____________________________________________________________________
701
702 Bool_t AliITSQADataMakerRec::AreEqual(Double_t a1,Double_t a2)
703 {
704   const Double_t kEpsilon= 1.e-14;
705   return TMath::Abs(a1-a2)<=kEpsilon*TMath::Abs(a1);      
706 }
707
708 //_________________________________________________________________
709 Int_t AliITSQADataMakerRec::GetDetTaskHisto(Int_t subdet,AliQAv1::TASKINDEX_t task)
710 {
711
712
713   Int_t histo=0;
714   switch(subdet)
715     {
716     case 1:
717       histo=fSPDDataMaker->GetTaskHisto(task);
718       //return histo;
719       break;
720     case 2:
721       histo=fSDDDataMaker->GetTaskHisto(task);
722       //return histo;
723       break;
724     case 3:
725       histo=fSSDDataMaker->GetTaskHisto(task);
726       //return histo;
727       break;
728     default:
729       AliWarning("No specific subdetector (SPD, SDD, SSD) selected!! Offset set to zero \n");
730       histo=0;
731       //return histo;
732       break;
733     }
734   //return offset;
735   return histo;
736 }