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