]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliGlobalQADataMaker.cxx
Re-organization of vertex constraints in the primary vertex determination:
[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(AliQA::TASKINDEX_t task, TObjArray * list)
22 {
23   //Detector specific actions at end of cycle
24   // do the QA checking
25   AliQAChecker::Instance()->Run(AliQA::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   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   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   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.02,6.30),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   Char_t *name[]={
146     "K0s mass (GeV)",
147     "Lambda0 + Lambda0Bar mass (GeV)"
148   };
149   Add2ESDsList(new TH1F(name[0],name[0],50,  0.4477,0.5477),kV0s0);
150   Add2ESDsList(new TH1F(name[1],name[1],50,  1.0657,1.1657),kV0s1);
151   }
152
153   {// PID related QA
154   Char_t *name[]={
155     "ITS: dEdx (ADC) for particles with momentum 0.4 - 0.5 (GeV)",
156     "TPC: dEdx (ADC) for particles with momentum 0.4 - 0.5 (GeV)",
157     "TOF: tracking - measured (ps)"
158   };
159   Add2ESDsList(new TH1F(name[0],name[0],50,0.00,200.),kPid0);
160   Add2ESDsList(new TH1F(name[1],name[1],50,0.00,100.),kPid1);
161   Add2ESDsList(new TH1F(name[2],name[2],50,-3500.,3500.),kPid2);
162   }
163
164 }
165
166 //____________________________________________________________________________
167 void AliGlobalQADataMaker::MakeRaws(AliRawReader* rawReader)
168 {
169   //Fill prepared histograms with Raw digit properties
170   rawReader->Reset() ;
171
172 }
173
174 //____________________________________________________________________________ 
175 void AliGlobalQADataMaker::MakeESDs(AliESDEvent * event) {
176   //-----------------------------------------------------------
177   // This function fills the ESD QA histograms
178   // as a part of global QA
179   //-----------------------------------------------------------
180   const AliESDEvent *esd=event;
181
182   Int_t ntrk=esd->GetNumberOfTracks() ; 
183   for (Int_t i=0; i<ntrk; i++) {
184     const AliESDtrack *track=esd->GetTrack(i);
185
186     // Cluster related QA
187     if (track->IsOn(AliESDtrack::kITSrefit)) {
188       Int_t n=track->GetITSclusters(0);
189       GetESDsData(kClr0)->Fill(Float_t(n)/6.); //6 is the number of ITS layers
190     }
191
192     if (track->IsOn(AliESDtrack::kTPCrefit)) {
193       Int_t n =track->GetTPCNcls();
194       Int_t nf=track->GetTPCNclsF();      // number of crossed TPC pad rows
195       if (nf>0) GetESDsData(kClr1)->Fill(Float_t(n)/nf);
196     }
197
198     if (track->IsOn(AliESDtrack::kTRDrefit)) {
199       Int_t n=track->GetTRDclusters(0);
200       GetESDsData(kClr2)->Fill(Float_t(n)/(6*24));//(6*24) is the number of TRD time bins
201     }
202
203     Double_t p=track->GetP();
204
205     // Track related QA
206     if (track->IsOn(AliESDtrack::kTPCrefit)) {
207       Float_t dz[2]; 
208       track->GetDZ(0.,0.,0.,esd->GetMagneticField(),dz); 
209       if ((TMath::Abs(dz[0])<3.) && (TMath::Abs(dz[1])<3.)) { // beam pipe
210         Double_t phi=track->Phi();
211         GetESDsData(kTrk0)->Fill(phi);
212         Double_t y=track->Eta();
213         GetESDsData(kTrk1)->Fill(y);
214
215         if (TMath::Abs(y)<0.9) {
216            GetESDsData(kTrk2)->Fill(p);
217            if (track->IsOn(AliESDtrack::kITSrefit)) GetESDsData(kTrk3)->Fill(p);
218            if (track->IsOn(AliESDtrack::kTOFout)) GetESDsData(kTrk4)->Fill(p);
219         }
220       }
221     }
222
223     // PID related QA
224     if ((p>0.4) && (p<0.5)) {
225       if (track->IsOn(AliESDtrack::kITSpid)) {
226         Double_t dedx=track->GetITSsignal();
227         GetESDsData(kPid0)->Fill(dedx);
228       }
229       if (track->IsOn(AliESDtrack::kTPCpid)) {
230         Double_t dedx=track->GetTPCsignal();
231         GetESDsData(kPid1)->Fill(dedx);
232       }
233     }
234     if (p>1.0) {
235       if (track->IsOn(AliESDtrack::kTOFpid)) {
236         Double_t times[10];
237         track->GetIntegratedTimes(times);
238         Double_t tof=track->GetTOFsignal();
239         GetESDsData(kPid2)->Fill(times[2]-tof);
240       }
241     }
242   }
243
244   TH1 *tpc=GetESDsData(kTrk2); tpc->Sumw2();
245   TH1 *its=GetESDsData(kTrk3); its->Sumw2();
246   TH1 *tof=GetESDsData(kTrk4); tof->Sumw2();
247   GetESDsData(kTrk5)->Divide(its,tpc,1,1.,"b");
248   GetESDsData(kTrk6)->Divide(tof,tpc,1,1.,"b");
249
250   // V0 related QA
251   Int_t nV0=esd->GetNumberOfV0s();
252   for (Int_t i=0; i<nV0; i++) {
253     Double_t mass;
254     AliESDv0 v0(*esd->GetV0(i));
255
256     v0.ChangeMassHypothesis(kK0Short);
257     mass=v0.GetEffMass();
258     GetESDsData(kV0s0)->Fill(mass);
259
260     v0.ChangeMassHypothesis(kLambda0);
261     mass=v0.GetEffMass();
262     GetESDsData(kV0s1)->Fill(mass);
263
264     v0.ChangeMassHypothesis(kLambda0Bar);
265     mass=v0.GetEffMass();
266     GetESDsData(kV0s1)->Fill(mass);
267   }
268
269 }