]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliGlobalQADataMaker.cxx
new classes for resonance and V0 analysis (from R. Vernet)
[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
62 void AliGlobalQADataMaker::InitESDs() {
63   //------------------------------------------------------
64   // This function books the ESD QA histograms
65   // as a part of global QA
66   //------------------------------------------------------
67
68   {// Cluster related QA
69   Char_t *name[]={
70     "Fraction of the assigned clusters in ITS",
71     "Fraction of the assigned clusters in TPC",
72     "Fraction of the assigned clusters in TRD"
73   };
74   Add2ESDsList(new TH1F(name[0],name[0],100,0.,2.),kClr0);
75   Add2ESDsList(new TH1F(name[1],name[1],100,0.,2.),kClr1);
76   Add2ESDsList(new TH1F(name[2],name[2],100,0.,2.),kClr2);
77   }
78
79   {// Track related QA
80   Char_t *name[]={
81     "Track azimuthal distribution (rad)",                   // kTrk0
82     "Track pseudo-rapidity distribution",                   // kTrk1
83     "TPC: track momentum distribution (GeV/c)",             // kTrk2
84     "TPC-ITS matched: track momentum distribution (GeV/c)", // kTrk3
85     "TPC-TOF matched: track momentum distribution (GeV/c)", // kTrk4
86     "TPC-ITS track-matching probability",                   // kTrk5
87     "TPC-TOF track-matching probability"                    // kTrk6
88   };
89   Add2ESDsList(new TH1F(name[0],name[0],100,-0.02,6.30),kTrk0);
90   Add2ESDsList(new TH1F(name[1],name[1],100,-2.00,2.00),kTrk1);
91   Add2ESDsList(new TH1F(name[2],name[2],50,  0.20,5.00),kTrk2);
92   Add2ESDsList(new TH1F(name[3],name[3],50,  0.20,5.00),kTrk3);
93   Add2ESDsList(new TH1F(name[4],name[4],50,  0.20,5.00),kTrk4);
94   Add2ESDsList(new TH1F(name[5],name[5],50,  0.20,5.00),kTrk5);
95   Add2ESDsList(new TH1F(name[6],name[6],50,  0.20,5.00),kTrk6);
96   }
97
98   {// V0 related QA
99   Char_t *name[]={
100     "K0s mass (GeV)",
101     "Lambda0 + Lambda0Bar mass (GeV)"
102   };
103   Add2ESDsList(new TH1F(name[0],name[0],50,  0.4477,0.5477),kV0s0);
104   Add2ESDsList(new TH1F(name[1],name[1],50,  1.0657,1.1657),kV0s1);
105   }
106
107   {// PID related QA
108   Char_t *name[]={
109     "ITS: dE/dx (ADC) for particles with momentum 0.4 - 0.5 (GeV/c)",
110     "TPC: dE/dx (ADC) for particles with momentum 0.4 - 0.5 (GeV/c)",
111     "TOF: tracking - measured (ps)"
112   };
113   Add2ESDsList(new TH1F(name[0],name[0],50,0.00,200.),kPid0);
114   Add2ESDsList(new TH1F(name[1],name[1],50,0.00,100.),kPid1);
115   Add2ESDsList(new TH1F(name[2],name[2],50,-3500.,3500.),kPid2);
116   }
117
118 }
119
120 void AliGlobalQADataMaker::MakeESDs(AliESDEvent * event) {
121   //-----------------------------------------------------------
122   // This function fills the ESD QA histograms
123   // as a part of global QA
124   //-----------------------------------------------------------
125   const AliESDEvent *esd=event;
126
127   Int_t ntrk=esd->GetNumberOfTracks() ; 
128   for (Int_t i=0; i<ntrk; i++) {
129     const AliESDtrack *track=esd->GetTrack(i);
130
131     // Cluster related QA
132     if (track->IsOn(AliESDtrack::kITSrefit)) {
133       Int_t n=track->GetITSclusters(0);
134       GetESDsData(kClr0)->Fill(Float_t(n)/6.); //6 is the number of ITS layers
135     }
136
137     if (track->IsOn(AliESDtrack::kTPCrefit)) {
138       Int_t n =track->GetTPCNcls();
139       Int_t nf=track->GetTPCNclsF();      // number of crossed TPC pad rows
140       GetESDsData(kClr1)->Fill(Float_t(n)/nf);
141     }
142
143     if (track->IsOn(AliESDtrack::kTRDrefit)) {
144       Int_t n=track->GetTRDclusters(0);
145       GetESDsData(kClr2)->Fill(Float_t(n)/120.);//120 is the number of TRD time bins
146     }
147
148     Double_t p=track->GetP();
149
150     // Track related QA
151     if (track->IsOn(AliESDtrack::kTPCrefit)) {
152       Float_t dz[2]; 
153       track->GetDZ(0.,0.,0.,esd->GetMagneticField(),dz); 
154       if ((TMath::Abs(dz[0])<3.) && (TMath::Abs(dz[1])<3.)) { // beam pipe
155         Double_t phi=track->Phi();
156         GetESDsData(kTrk0)->Fill(phi);
157         Double_t y=track->Y();
158         GetESDsData(kTrk1)->Fill(y);
159
160         GetESDsData(kTrk2)->Fill(p);
161
162         if (track->IsOn(AliESDtrack::kITSrefit)) 
163            GetESDsData(kTrk3)->Fill(p);
164
165         if (TMath::Abs(y)<0.9)
166            if (track->IsOn(AliESDtrack::kTOFout))
167               GetESDsData(kTrk4)->Fill(p);
168       }
169     }
170
171     // PID related QA
172     if ((p>0.4) && (p<0.5)) {
173       if (track->IsOn(AliESDtrack::kITSpid)) {
174         Double_t dedx=track->GetITSsignal();
175         GetESDsData(kPid0)->Fill(dedx);
176       }
177       if (track->IsOn(AliESDtrack::kTPCpid)) {
178         Double_t dedx=track->GetTPCsignal();
179         GetESDsData(kPid1)->Fill(dedx);
180       }
181     }
182     if (p>1.0) {
183       if (track->IsOn(AliESDtrack::kTOFpid)) {
184         Double_t times[10];
185         track->GetIntegratedTimes(times);
186         Double_t tof=track->GetTOFsignal();
187         GetESDsData(kPid2)->Fill(times[2]-tof);
188       }
189     }
190   }
191
192   TH1 *tpc=GetESDsData(kTrk2);
193   TH1 *its=GetESDsData(kTrk3);
194   TH1 *tof=GetESDsData(kTrk4);
195   GetESDsData(kTrk5)->Divide(its,tpc,1,1.,"b");
196   GetESDsData(kTrk6)->Divide(tof,tpc,1,1.,"b");
197
198   // V0 related QA
199   Int_t nV0=esd->GetNumberOfV0s();
200   for (Int_t i=0; i<nV0; i++) {
201     Double_t mass;
202     AliESDv0 v0(*esd->GetV0(i));
203
204     v0.ChangeMassHypothesis(kK0Short);
205     mass=v0.GetEffMass();
206     GetESDsData(kV0s0)->Fill(mass);
207
208     v0.ChangeMassHypothesis(kLambda0);
209     mass=v0.GetEffMass();
210     GetESDsData(kV0s1)->Fill(mass);
211
212     v0.ChangeMassHypothesis(kLambda0Bar);
213     mass=v0.GetEffMass();
214     GetESDsData(kV0s1)->Fill(mass);
215   }
216
217 }