]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliGlobalQADataMaker.cxx
Bug fix
[u/mrichter/AliRoot.git] / STEER / AliGlobalQADataMaker.cxx
1 /*
2  The class for calculating the global (not detector specific) quality assurance.
3  It reuses the following TLists from its base class 
4     AliQADataMaker::fRecPointsQAList (for keeping the track residuals)
5     AliQADataMaker::fESDsQAList      (for keeping global ESD QA data)
6 */
7
8 #include <TH1F.h>
9
10 #include "AliGlobalQADataMaker.h"
11 #include "AliGeomManager.h"
12 #include "AliESDEvent.h"
13
14 ClassImp(AliGlobalQADataMaker)
15  
16 void AliGlobalQADataMaker::InitRecPoints() {
17   //------------------------------------------------------
18   // This function books the histograms of *track*residuals*
19   // as a part of global QA
20   //------------------------------------------------------
21   Char_t *name[]={
22     "SPD1 residuals Y","SPD1 residuals Z",
23     "SPD2 residuals Y","SPD2 residuals Z",
24     "SDD1 residuals Y","SDD1 residuals Z",
25     "SDD2 residuals Y","SDD2 residuals Z",
26     "SSD1 residuals Y","SSD1 residuals Z",
27     "SSD2 residuals Y","SSD2 residuals Z",
28
29     "TPC1 residuals Y","TPC1 residuals Z",
30     "TPC2 residuals Y","TPC2 residuals Z",
31
32     "TRD1 residuals Y","TRD1 residuals Z",
33     "TRD2 residuals Y","TRD2 residuals Z",
34     "TRD3 residuals Y","TRD3 residuals Z",
35     "TRD4 residuals Y","TRD4 residuals Z",
36     "TRD5 residuals Y","TRD5 residuals Z",
37     "TRD6 residuals Y","TRD6 residuals Z",
38
39     "TOF residuals Y","TOF residuals Z",
40
41     "PHOS1 residuals Y","PHOS1 residuals Z",
42     "PHOS2 residuals Y","PHOS2 residuals Z",
43
44     "HMPID residuals Y","HMPID residuals Z",
45
46     "MUON residuals Y","MUON residuals Z",
47
48     "EMCAL residuals Y","EMCAL residuals Z"
49   };
50
51   for (Int_t m=1; m<AliGeomManager::kLastLayer; m++) {
52     Int_t i=2*m-2;
53     TH1F *h=new TH1F(name[i],name[i],100,-5.,5.);
54     Add2RecPointsList(h,i);    
55     h=new TH1F(name[i+1],name[i+1],100,-5.,5.);
56     Add2RecPointsList(h,i+1);    
57   }
58 }
59
60 void AliGlobalQADataMaker::InitESDs() {
61   //------------------------------------------------------
62   // This function books the ESD QA histograms
63   // as a part of global QA
64   //------------------------------------------------------
65   Char_t *name[]={
66     "Fraction of the assigned clusters in ITS",
67     "Fraction of the assigned clusters in TPC",
68     "Fraction of the assigned clusters in TRD"
69   };
70   Add2ESDsList(new TH1F(name[0],name[0],100,0.,2.),0);
71   Add2ESDsList(new TH1F(name[1],name[1],100,0.,2.),1);
72   Add2ESDsList(new TH1F(name[2],name[2],100,0.,2.),2);
73 }
74
75 void AliGlobalQADataMaker::MakeESDs(AliESDEvent * esd) {
76   //-----------------------------------------------------------
77   // This function fills the ESD QA histograms
78   // as a part of global QA
79   //-----------------------------------------------------------
80   Int_t ntrk=esd->GetNumberOfTracks() ; 
81   for (Int_t i=0; i<ntrk; i++) {
82     AliESDtrack *track=esd->GetTrack(i);
83
84     if (track->IsOn(AliESDtrack::kITSrefit)) {
85       Int_t n=track->GetITSclusters(0);
86       GetESDsData(0)->Fill(Float_t(n)/6.); //6 is the number of ITS layers
87     }
88
89     if (track->IsOn(AliESDtrack::kTPCrefit)) {
90       Int_t n =track->GetTPCNcls();
91       Int_t nf=track->GetTPCNclsF();      // number of crossed TPC pad rows
92       GetESDsData(1)->Fill(Float_t(n)/nf);
93     }
94
95     if (track->IsOn(AliESDtrack::kTRDrefit)) {
96       Int_t n=track->GetTRDclusters(0);
97       GetESDsData(2)->Fill(Float_t(n)/120.);//120 is the number of TRD time bins
98     }
99
100   }
101 }