]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliGlobalQADataMaker.cxx
Adding also the SSD absolute residuals in Y
[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 <TPDGCode.h>
9 #include <TH1F.h>
10
11 #include "AliGlobalQADataMaker.h"
12 #include "AliGeomManager.h"
13 #include "AliESDEvent.h"
14 #include "AliESDv0.h"
15
16 ClassImp(AliGlobalQADataMaker)
17  
18 void AliGlobalQADataMaker::InitRecPoints() {
19   //------------------------------------------------------
20   // This function books the histograms of *track*residuals*
21   // as a part of global QA
22   //------------------------------------------------------
23   Char_t *name[]={
24     "SPD1 residuals Y","SPD1 residuals Z",
25     "SPD2 residuals Y","SPD2 residuals Z",
26     "SDD1 residuals Y","SDD1 residuals Z",
27     "SDD2 residuals Y","SDD2 residuals Z",
28     "SSD1 residuals Y","SSD1 residuals Z",
29     "SSD2 residuals Y","SSD2 residuals Z",
30
31     "TPC1 residuals Y","TPC1 residuals Z",
32     "TPC2 residuals Y","TPC2 residuals Z",
33
34     "TRD1 residuals Y","TRD1 residuals Z",
35     "TRD2 residuals Y","TRD2 residuals Z",
36     "TRD3 residuals Y","TRD3 residuals Z",
37     "TRD4 residuals Y","TRD4 residuals Z",
38     "TRD5 residuals Y","TRD5 residuals Z",
39     "TRD6 residuals Y","TRD6 residuals Z",
40
41     "TOF residuals Y","TOF residuals Z",
42
43     "PHOS1 residuals Y","PHOS1 residuals Z",
44     "PHOS2 residuals Y","PHOS2 residuals Z",
45
46     "HMPID residuals Y","HMPID residuals Z",
47
48     "MUON residuals Y","MUON residuals Z",
49
50     "EMCAL residuals Y","EMCAL residuals Z"
51   };
52
53   for (Int_t m=1; m<AliGeomManager::kLastLayer; m++) {
54     Int_t i=2*m-2;
55     TH1F *h=new TH1F(name[i],name[i],100,-5.,5.);
56     Add2RecPointsList(h,i);    
57     h=new TH1F(name[i+1],name[i+1],100,-5.,5.);
58     Add2RecPointsList(h,i+1);    
59   }
60
61   Add2RecPointsList(
62   new TH1F("SSD1 absolute residuals Y for Z<0 (cm)",
63            "SSD1 absolute residuals Y for Z<0 (cm)",100,-2.,2.),40);
64   Add2RecPointsList(
65   new TH1F("SSD1 absolute residuals Z for Z<0 (cm)",
66            "SSD1 absolute residuals Z for Z<0 (cm)",100,-2.,2.),41);
67   Add2RecPointsList(
68   new TH1F("SSD1 absolute residuals Y for Z>0 (cm)",
69            "SSD1 absolute residuals Y for Z>0 (cm)",100,-2.,2.),42);
70   Add2RecPointsList(
71   new TH1F("SSD1 absolute residuals Z for Z>0 (cm)",
72            "SSD1 absolute residuals Z for Z>0 (cm)",100,-2.,2.),43);
73
74
75   Add2RecPointsList(
76   new TH1F("SSD2 absolute residuals Y for Z<0 (cm)",
77            "SSD2 absolute residuals Y for Z<0 (cm)",100,-3.,3.),44);
78   Add2RecPointsList(
79   new TH1F("SSD2 absolute residuals Z for Z<0 (cm)",
80            "SSD2 absolute residuals Z for Z<0 (cm)",100,-3.,3.),45);
81   Add2RecPointsList(
82   new TH1F("SSD2 absolute residuals Y for Z>0 (cm)",
83            "SSD2 absolute residuals Y for Z>0 (cm)",100,-3.,3.),46);
84   Add2RecPointsList(
85   new TH1F("SSD2 absolute residuals Z for Z>0 (cm)",
86            "SSD2 absolute residuals Z for Z>0 (cm)",100,-3.,3.),47);
87   
88 }
89
90 void AliGlobalQADataMaker::InitESDs() {
91   //------------------------------------------------------
92   // This function books the ESD QA histograms
93   // as a part of global QA
94   //------------------------------------------------------
95
96   {// Cluster related QA
97   Char_t *name[]={
98     "Fraction of the assigned clusters in ITS",
99     "Fraction of the assigned clusters in TPC",
100     "Fraction of the assigned clusters in TRD"
101   };
102   Add2ESDsList(new TH1F(name[0],name[0],100,0.,2.),kClr0);
103   Add2ESDsList(new TH1F(name[1],name[1],100,0.,2.),kClr1);
104   Add2ESDsList(new TH1F(name[2],name[2],100,0.,2.),kClr2);
105   }
106
107   {// Track related QA
108   Char_t *name[]={
109     "Track azimuthal distribution (rad)",                   // kTrk0
110     "Track pseudo-rapidity distribution",                   // kTrk1
111     "TPC: track momentum distribution (GeV)",               // kTrk2
112     "TPC-ITS matched: track momentum distribution (GeV)",   // kTrk3
113     "TPC-TOF matched: track momentum distribution (GeV)",   // kTrk4
114     "TPC-ITS track-matching probability",                   // kTrk5
115     "TPC-TOF track-matching probability"                    // kTrk6
116   };
117   Add2ESDsList(new TH1F(name[0],name[0],100,-0.02,6.30),kTrk0);
118   Add2ESDsList(new TH1F(name[1],name[1],100,-2.00,2.00),kTrk1);
119   Add2ESDsList(new TH1F(name[2],name[2],50,  0.20,5.00),kTrk2);
120   Add2ESDsList(new TH1F(name[3],name[3],50,  0.20,5.00),kTrk3);
121   Add2ESDsList(new TH1F(name[4],name[4],50,  0.20,5.00),kTrk4);
122   Add2ESDsList(new TH1F(name[5],name[5],50,  0.20,5.00),kTrk5);
123   Add2ESDsList(new TH1F(name[6],name[6],50,  0.20,5.00),kTrk6);
124   }
125
126   {// V0 related QA
127   Char_t *name[]={
128     "K0s mass (GeV)",
129     "Lambda0 + Lambda0Bar mass (GeV)"
130   };
131   Add2ESDsList(new TH1F(name[0],name[0],50,  0.4477,0.5477),kV0s0);
132   Add2ESDsList(new TH1F(name[1],name[1],50,  1.0657,1.1657),kV0s1);
133   }
134
135   {// PID related QA
136   Char_t *name[]={
137     "ITS: dEdx (ADC) for particles with momentum 0.4 - 0.5 (GeV)",
138     "TPC: dEdx (ADC) for particles with momentum 0.4 - 0.5 (GeV)",
139     "TOF: tracking - measured (ps)"
140   };
141   Add2ESDsList(new TH1F(name[0],name[0],50,0.00,200.),kPid0);
142   Add2ESDsList(new TH1F(name[1],name[1],50,0.00,100.),kPid1);
143   Add2ESDsList(new TH1F(name[2],name[2],50,-3500.,3500.),kPid2);
144   }
145
146 }
147
148 void AliGlobalQADataMaker::MakeESDs(AliESDEvent * event) {
149   //-----------------------------------------------------------
150   // This function fills the ESD QA histograms
151   // as a part of global QA
152   //-----------------------------------------------------------
153   const AliESDEvent *esd=event;
154
155   Int_t ntrk=esd->GetNumberOfTracks() ; 
156   for (Int_t i=0; i<ntrk; i++) {
157     const AliESDtrack *track=esd->GetTrack(i);
158
159     // Cluster related QA
160     if (track->IsOn(AliESDtrack::kITSrefit)) {
161       Int_t n=track->GetITSclusters(0);
162       GetESDsData(kClr0)->Fill(Float_t(n)/6.); //6 is the number of ITS layers
163     }
164
165     if (track->IsOn(AliESDtrack::kTPCrefit)) {
166       Int_t n =track->GetTPCNcls();
167       Int_t nf=track->GetTPCNclsF();      // number of crossed TPC pad rows
168       GetESDsData(kClr1)->Fill(Float_t(n)/nf);
169     }
170
171     if (track->IsOn(AliESDtrack::kTRDrefit)) {
172       Int_t n=track->GetTRDclusters(0);
173       GetESDsData(kClr2)->Fill(Float_t(n)/(6*24));//(6*24) is the number of TRD time bins
174     }
175
176     Double_t p=track->GetP();
177
178     // Track related QA
179     if (track->IsOn(AliESDtrack::kTPCrefit)) {
180       Float_t dz[2]; 
181       track->GetDZ(0.,0.,0.,esd->GetMagneticField(),dz); 
182       if ((TMath::Abs(dz[0])<3.) && (TMath::Abs(dz[1])<3.)) { // beam pipe
183         Double_t phi=track->Phi();
184         GetESDsData(kTrk0)->Fill(phi);
185         Double_t y=track->Eta();
186         GetESDsData(kTrk1)->Fill(y);
187
188         if (TMath::Abs(y)<0.9) {
189            GetESDsData(kTrk2)->Fill(p);
190            if (track->IsOn(AliESDtrack::kITSrefit)) GetESDsData(kTrk3)->Fill(p);
191            if (track->IsOn(AliESDtrack::kTOFout)) GetESDsData(kTrk4)->Fill(p);
192         }
193       }
194     }
195
196     // PID related QA
197     if ((p>0.4) && (p<0.5)) {
198       if (track->IsOn(AliESDtrack::kITSpid)) {
199         Double_t dedx=track->GetITSsignal();
200         GetESDsData(kPid0)->Fill(dedx);
201       }
202       if (track->IsOn(AliESDtrack::kTPCpid)) {
203         Double_t dedx=track->GetTPCsignal();
204         GetESDsData(kPid1)->Fill(dedx);
205       }
206     }
207     if (p>1.0) {
208       if (track->IsOn(AliESDtrack::kTOFpid)) {
209         Double_t times[10];
210         track->GetIntegratedTimes(times);
211         Double_t tof=track->GetTOFsignal();
212         GetESDsData(kPid2)->Fill(times[2]-tof);
213       }
214     }
215   }
216
217   TH1 *tpc=GetESDsData(kTrk2); tpc->Sumw2();
218   TH1 *its=GetESDsData(kTrk3); its->Sumw2();
219   TH1 *tof=GetESDsData(kTrk4); tof->Sumw2();
220   GetESDsData(kTrk5)->Divide(its,tpc,1,1.,"b");
221   GetESDsData(kTrk6)->Divide(tof,tpc,1,1.,"b");
222
223   // V0 related QA
224   Int_t nV0=esd->GetNumberOfV0s();
225   for (Int_t i=0; i<nV0; i++) {
226     Double_t mass;
227     AliESDv0 v0(*esd->GetV0(i));
228
229     v0.ChangeMassHypothesis(kK0Short);
230     mass=v0.GetEffMass();
231     GetESDsData(kV0s0)->Fill(mass);
232
233     v0.ChangeMassHypothesis(kLambda0);
234     mass=v0.GetEffMass();
235     GetESDsData(kV0s1)->Fill(mass);
236
237     v0.ChangeMassHypothesis(kLambda0Bar);
238     mass=v0.GetEffMass();
239     GetESDsData(kV0s1)->Fill(mass);
240   }
241
242 }