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