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