]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliGlobalQADataMaker.cxx
A better binning for the new histogramms.
[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 #include <TH2F.h>
11
12 #include "AliQAChecker.h"
13 #include "AliGlobalQADataMaker.h"
14 #include "AliGeomManager.h"
15 #include "AliESDEvent.h"
16 #include "AliESDv0.h"
17 #include "AliRawReader.h"
18 #include "AliESDVZERO.h"
19 #include "AliMultiplicity.h" 
20
21 ClassImp(AliGlobalQADataMaker)
22  
23 //____________________________________________________________________________ 
24 void AliGlobalQADataMaker::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
25 {
26   //Detector specific actions at end of cycle
27   // do the QA checking
28   AliQAChecker::Instance()->Run(AliQAv1::kGLOBAL, task, list) ;  
29 }
30
31 //____________________________________________________________________________ 
32 void AliGlobalQADataMaker::InitRaws()
33 {
34   // create Raws histograms in Raws subdir
35 }
36
37 //____________________________________________________________________________ 
38 void AliGlobalQADataMaker::InitRecPoints() {
39   //------------------------------------------------------
40   // This function books the histograms of *track*residuals*
41   // as a part of global QA
42   //------------------------------------------------------
43   const Char_t *name[]={
44     "hGlobalSPD1ResidualsY","SPD1ResidualsZ",
45     "hGlobalSPD2ResidualsY","SPD2ResidualsZ",
46     "hGlobalSDD1ResidualsY","SDD1ResidualsZ",
47     "hGlobalSDD2ResidualsY","SDD2ResidualsZ",
48     "hGlobalSSD1ResidualsY","SSD1ResidualsZ",
49     "hGlobalSSD2ResidualsY","SSD2ResidualsZ",
50     
51     "hGlobalTPC1ResidualsY","TPC1ResidualsZ",
52     "hGlobalTPC2ResidualsY","TPC2ResidualsZ",
53     
54     "hGlobalTRD1ResidualsY","TRD1ResidualsZ",
55     "hGlobalTRD2ResidualsY","TRD2ResidualsZ",
56     "hGlobalTRD3ResidualsY","TRD3ResidualsZ",
57     "hGlobalTRD4ResidualsY","TRD4ResidualsZ",
58     "hGlobalTRD5ResidualsY","TRD5ResidualsZ",
59     "hGlobalTRD6ResidualsY","TRD6ResidualsZ",
60     
61     "hGlobalTOFResidualsY","TOFResidualsZ",
62     
63     "hGlobalPHOS1ResidualsY","PHOS1ResidualsZ",
64     "hGlobalPHOS2ResidualsY","PHOS2ResidualsZ",
65     
66     "hGlobalHMPIDResidualsY","HMPIDResidualsZ",
67     
68     "hGlobalMUONResidualsY","MUONResidualsZ",
69     
70     "hGlobalEMCALResidualsY","EMCALResidualsZ"
71   };
72   const Char_t *title[]={
73     "SPD1 residuals Y","SPD1 residuals Z",
74     "SPD2 residuals Y","SPD2 residuals Z",
75     "SDD1 residuals Y","SDD1 residuals Z",
76     "SDD2 residuals Y","SDD2 residuals Z",
77     "SSD1 residuals Y","SSD1 residuals Z",
78     "SSD2 residuals Y","SSD2 residuals Z",
79     
80     "TPC1 residuals Y","TPC1 residuals Z",
81     "TPC2 residuals Y","TPC2 residuals Z",
82     
83     "TRD1 residuals Y","TRD1 residuals Z",
84     "TRD2 residuals Y","TRD2 residuals Z",
85     "TRD3 residuals Y","TRD3 residuals Z",
86     "TRD4 residuals Y","TRD4 residuals Z",
87     "TRD5 residuals Y","TRD5 residuals Z",
88     "TRD6 residuals Y","TRD6 residuals Z",
89     
90     "TOF residuals Y","TOF residuals Z",
91     
92     "PHOS1 residuals Y","PHOS1 residuals Z",
93     "PHOS2 residuals Y","PHOS2 residuals Z",
94     
95     "HMPID residuals Y","HMPID residuals Z",
96     
97     "MUON residuals Y","MUON residuals Z",
98     
99     "EMCAL residuals Y","EMCAL residuals Z"
100   };
101   
102   for (Int_t m=1; m<AliGeomManager::kLastLayer; m++) {
103     Int_t i=2*m-2;
104     TH1F *h=new TH1F(name[i],title[i],100,-5.,5.);
105     Add2RecPointsList(h,i);    
106     h=new TH1F(name[i+1],title[i+1],100,-5.,5.);
107     Add2RecPointsList(h,i+1);    
108   }
109
110   Add2RecPointsList(
111   new TH1F("hGlobalSSD1AbsoluteResidualsYNegZ",
112            "SSD1 Absolute Residuals Y Neg Z",100,-2.,2.),40);
113   Add2RecPointsList(
114   new TH1F("hGlobalSSD1AbsoluteResidualsZNegZ",
115            "SSD1 Absolute Residuals Z Neg Z",100,-2.,2.),41);
116   Add2RecPointsList(
117   new TH1F("hGlobalSSD1AbsoluteResidualsYPosZ",
118            "SSD1 Absolute Residuals Y Pos Z",100,-2.,2.),42);
119   Add2RecPointsList(
120   new TH1F("hGlobalSSD1AbsoluteResidualsZPosZ",
121            "SSD1 Absolute Residuals Z Pos Z",100,-2.,2.),43);
122
123
124   Add2RecPointsList(
125   new TH1F("hGlobalSSD2AbsoluteResidualsYNegZ",
126            "SSD2 Absolute Residuals Y Neg Z",100,-3.,3.),44);
127   Add2RecPointsList(
128   new TH1F("hGlobalSSD2AbsoluteResidualsZNegZ",
129            "SSD2 Absolute Residuals Z Neg Z",100,-3.,3.),45);
130   Add2RecPointsList(
131   new TH1F("hGlobalSSD2AbsoluteResidualsYPosZ",
132            "SSD2 Absolute Residuals Y Pos Z",100,-3.,3.),46);
133   Add2RecPointsList(
134   new TH1F("hGlobalSSD2AbsoluteResidualsZPosZ",
135            "SSD2Absolute Residuals Z Pos Z",100,-3.,3.),47);
136   
137 }
138
139 //____________________________________________________________________________ 
140 void AliGlobalQADataMaker::InitESDs() {
141   //------------------------------------------------------
142   // This function books the ESD QA histograms
143   // as a part of global QA
144   //------------------------------------------------------
145
146   const Bool_t expert   = kTRUE ; 
147   const Bool_t image    = kTRUE ; 
148  {// Cluster related QA
149     const Char_t *name[]={
150       "hGlobalFractionAssignedClustersITS",
151       "hGlobalFractionAssignedClustersTPC",
152       "hGlobalFractionAssignedClustersTRD",
153       "hGlobalClustersPerITSModule"
154     };
155     const Char_t *title[]={
156       "Fraction of the assigned clusters in ITS",
157       "Fraction of the assigned clusters in TPC",
158       "Fraction of the assigned clusters in TRD",
159       "Number of clusters per an ITS module"
160     };
161     Add2ESDsList(new TH1F(name[0],title[0],100,0.,2.),kClr0, !expert, image);
162     Add2ESDsList(new TH1F(name[1],title[1],100,0.,2.),kClr1, !expert, image);
163     Add2ESDsList(new TH1F(name[2],title[2],100,0.,2.),kClr2, !expert, image);
164     Add2ESDsList(new TH1F(name[3],title[3],2201,-0.5,2200.5),kClr3, !expert, image);
165   }
166
167   {// Track related QA
168     const Char_t *name[]={
169       "hGlobalTrackAzimuthe",                               // kTrk0
170       "hGlobalTrackEta",                                    // kTrk1
171       "hGlobalTPCTrackpT",                                  // kTrk2
172       "hGlobalTPCITSMatchedpT",                             // kTrk3
173       "hGlobalTPCTOFMatchedpT",                             // kTrk4
174       "hGlobalTPCITSMatchingProbability",                   // kTrk5
175       "hGlobalTPCTOFMatchingProbability",                   // kTrk6
176       "hGlobalTPCsideAposDCA",                              // kTrk7
177       "hGlobalTPCsideAnegDCA",                              // kTrk8
178       "hGlobalTPCsideCposDCA",                              // kTrk9
179       "hGlobalTPCsideCnegDCA"                               // kTrk10
180   };
181     const Char_t *title[]={
182       "Track azimuthal distribution (rad)",                   // kTrk0
183       "Track pseudo-rapidity distribution",                   // kTrk1
184       "TPC: track momentum distribution (GeV)",               // kTrk2
185       "TPC-ITS matched: track momentum distribution (GeV)",   // kTrk3
186       "TPC-TOF matched: track momentum distribution (GeV)",   // kTrk4
187       "TPC-ITS track-matching probability",                   // kTrk5
188       "TPC-TOF track-matching probability",                   // kTrk6
189       "TPC side A: DCA for the positive tracks (mm)",         // kTrk7
190       "TPC side A: DCA for the negative tracks (mm)",         // kTrk8
191       "TPC side C: DCA for the positive tracks (mm)",         // kTrk9
192       "TPC side C: DCA for the negative tracks (mm)"          // kTrk10
193     };
194   Add2ESDsList(new TH1F(name[0],title[0],100, 0.,TMath::TwoPi()),kTrk0, !expert, image);
195   Add2ESDsList(new TH1F(name[1],title[1],100,-2.00,2.00),kTrk1, !expert, image);
196   Add2ESDsList(new TH1F(name[2],title[2],50,  0.20,5.00),kTrk2, !expert, image);
197   Add2ESDsList(new TH1F(name[3],title[3],50,  0.20,5.00),kTrk3, !expert, image);
198   Add2ESDsList(new TH1F(name[4],title[4],50,  0.20,5.00),kTrk4, !expert, image);
199   Add2ESDsList(new TH1F(name[5],title[5],50,  0.20,5.00),kTrk5, !expert, image);
200   Add2ESDsList(new TH1F(name[6],title[6],50,  0.20,5.00),kTrk6, !expert, image);
201   Add2ESDsList(new TH1F(name[7],title[7],50, -25.0,25.0),kTrk7, !expert, image);
202   Add2ESDsList(new TH1F(name[8],title[8],50, -25.0,25.0),kTrk8, !expert, image);
203   Add2ESDsList(new TH1F(name[9],title[9],50, -25.0,25.0),kTrk9, !expert, image);
204   Add2ESDsList(new TH1F(name[10],title[10],50, -25.0,25.0),kTrk10, !expert, image);
205   }
206
207   {// V0 related QA
208     const Char_t *name[]={
209       "hGlobalPromptK0sMass",
210       "hGlobalOfflineK0sMass",
211       "hGlobalPromptLambda0Lambda0BarMass",
212       "hGlobalOfflineLambda0Lambda0BarMass"
213     };
214     const Char_t *title[]={
215       "On-the-fly K0s mass (GeV)",
216       "Offline K0s mass (GeV)",
217       "On-the-fly Lambda0 + Lambda0Bar mass (GeV)",
218       "Offline Lambda0 + Lambda0Bar mass (GeV)"
219     };
220     Add2ESDsList(new TH1F(name[0],title[0],50,  0.4477,0.5477),kK0on, !expert, image);
221     Add2ESDsList(new TH1F(name[1],title[1],50,  0.4477,0.5477),kK0off, !expert, image);
222     Add2ESDsList(new TH1F(name[2],title[2],50,  1.0657,1.1657),kL0on, !expert, image);
223     Add2ESDsList(new TH1F(name[3],title[3],50,  1.0657,1.1657),kL0off, !expert, image);
224   }
225
226   {// PID related QA
227   const Char_t *name[]={
228     "hGlobalITSdEdx",
229     "hGlobalTPCdEdx",
230     "hGlobalTOFTrackingvsMeasured",
231     "hGlobalTPCdEdxvsMomentum"
232    };
233     const Char_t *title[]={
234       "ITS: dEdx (ADC) for particles with momentum 0.4 - 0.5 (GeV)",
235       "TPC: dEdx (ADC) for particles with momentum 0.4 - 0.5 (GeV)",
236       "TOF: tracking - measured (ps)",
237       "TPC: dEdx (A.U.) vs momentum (GeV)"
238      };
239     Add2ESDsList(new TH1F(name[0],title[0],50,0.00,200.),kPid0, !expert, image);
240     Add2ESDsList(new TH1F(name[1],title[1],50,0.00,100.),kPid1, !expert, image);
241     Add2ESDsList(new TH1F(name[2],title[2],50,-3500.,3500.),kPid2, !expert, image);
242     Add2ESDsList(new TH2F(name[3],title[3],1500,0.05,15.,700,0.,700.),kPid3,!expert,image);
243    }
244   {// Multiplicity related QA
245     const Char_t *name[]={
246       "hGlobalV0AvsITS",
247       "hGlobalV0CvsITS"
248     };
249     const Char_t *title[]={
250       "Multiplicity: V0A vs ITS",
251       "Multiplicity: V0C vs ITS"
252     };
253     TH2F *h0=new TH2F(name[0],title[0],41,-0.5,40.5, 33,-0.5,32.5);
254     Add2ESDsList(h0,kMlt0, !expert, image);
255     TH2F *h1=new TH2F(name[1],title[1],41,-0.5,40.5, 33,-0.5,32.5);
256     Add2ESDsList(h1,kMlt1, !expert, image);
257   }
258
259 }
260
261 //____________________________________________________________________________
262 void AliGlobalQADataMaker::MakeRaws(AliRawReader* rawReader)
263 {
264   //Fill prepared histograms with Raw digit properties
265   rawReader->Reset() ;
266
267 }
268
269 //____________________________________________________________________________ 
270 void AliGlobalQADataMaker::MakeESDs(AliESDEvent * event) {
271   //-----------------------------------------------------------
272   // This function fills the ESD QA histograms
273   // as a part of global QA
274   //-----------------------------------------------------------
275   // Check id histograms already created for this Event Specie
276   if ( ! GetESDsData(kClr0) )
277     InitESDs() ;
278
279
280   const AliESDEvent *esd=event;
281
282   Int_t ntrk=esd->GetNumberOfTracks() ; 
283   for (Int_t i=0; i<ntrk; i++) {
284     const AliESDtrack *track=esd->GetTrack(i);
285
286     // Cluster related QA
287     if (track->IsOn(AliESDtrack::kITSrefit)) {
288       Int_t n=track->GetITSclusters(0);
289       GetESDsData(kClr0)->Fill(Float_t(n)/6.); //6 is the number of ITS layers
290     }
291
292     for (Int_t i=0; i<6; i++) {
293       Int_t idet, sts;
294       Float_t xloc,zloc;
295       if (!track->GetITSModuleIndexInfo(i,idet,sts,xloc,zloc)) continue;
296       if (i>=2) idet+=240;
297       if (i>=4) idet+=260;
298       if ((sts==1)||(sts==2)||(sts==4)) GetESDsData(kClr3)->Fill(idet);  
299     }
300
301     if (track->IsOn(AliESDtrack::kTPCrefit)) {
302       Int_t n =track->GetTPCNcls();
303       Int_t nf=track->GetTPCNclsF();      // number of crossed TPC pad rows
304       if (nf>0) {
305         Double_t val = n*1.0/nf; 
306         GetESDsData(kClr1)->Fill(val); 
307       }
308     }
309
310     if (track->IsOn(AliESDtrack::kTRDrefit)) {
311       Int_t n=track->GetTRDclusters(0);
312       GetESDsData(kClr2)->Fill(Float_t(n)/(6*24));//(6*24) is the number of TRD time bins
313     }
314
315     Double_t p=track->GetP();
316
317     // Track related QA
318     if (track->IsOn(AliESDtrack::kTPCrefit)) {
319       Float_t dz[2]; 
320       track->GetDZ(0.,0.,0.,esd->GetMagneticField(),dz); 
321       if ((TMath::Abs(dz[0])<3.) && (TMath::Abs(dz[1])<3.)) { // beam pipe
322         Double_t phi=track->Phi();
323         GetESDsData(kTrk0)->Fill(phi);
324         Double_t y=track->Eta();
325         GetESDsData(kTrk1)->Fill(y);
326
327         if (TMath::Abs(y)<0.9) {
328            GetESDsData(kTrk2)->Fill(p);
329            if (track->IsOn(AliESDtrack::kITSrefit)) GetESDsData(kTrk3)->Fill(p);
330           //if (track->IsOn(AliESDtrack::kTOFout)) GetESDsData(kTrk4)->Fill(p);
331            if (track->GetTOFsignal()>0) GetESDsData(kTrk4)->Fill(p);
332         }
333       }
334     }
335     const AliExternalTrackParam *tpcTrack=track->GetTPCInnerParam();
336     const AliExternalTrackParam *innTrack=track->GetInnerParam();
337     if (tpcTrack)
338     if (innTrack) {
339        const AliESDVertex *vtx=esd->GetPrimaryVertex();
340        Double_t xv=vtx->GetXv();
341        Double_t yv=vtx->GetYv();
342        Double_t zv=vtx->GetZv();
343        Float_t dz[2];
344        tpcTrack->GetDZ(xv,yv,zv,esd->GetMagneticField(),dz);
345        dz[0]*=10.; // in mm
346        if (innTrack->GetZ()  > 0)
347        if (innTrack->GetTgl()> 0) { // TPC side A
348           if (tpcTrack->GetSign() > 0) GetESDsData(kTrk7)->Fill(dz[0]);
349           else                         GetESDsData(kTrk8)->Fill(dz[0]);
350        }
351        if (innTrack->GetZ()  < 0)
352        if (innTrack->GetTgl()< 0) { // TPC side C
353           if (tpcTrack->GetSign() > 0) GetESDsData(kTrk9)->Fill(dz[0]);
354           else                         GetESDsData(kTrk10)->Fill(dz[0]);
355        }
356     }
357
358     // PID related QA
359     if ((p>0.4) && (p<0.5)) {
360       if (track->IsOn(AliESDtrack::kITSpid)) {
361         Double_t dedx=track->GetITSsignal();
362         GetESDsData(kPid0)->Fill(dedx);
363       }
364       if (track->IsOn(AliESDtrack::kTPCpid)) {
365         Double_t dedx=track->GetTPCsignal();
366         GetESDsData(kPid1)->Fill(dedx);
367       }
368     }
369     if (p>1.0) {
370       if (track->IsOn(AliESDtrack::kTOFpid)) {
371         Double_t times[10];
372         track->GetIntegratedTimes(times);
373         Double_t tof=track->GetTOFsignal();
374         GetESDsData(kPid2)->Fill(times[2]-tof);
375       }
376     }
377     const AliExternalTrackParam *par=track->GetInnerParam();
378     if (par) {
379       Double_t pp=par->GetP();
380       Double_t dedx=track->GetTPCsignal();
381       TH2F *h = dynamic_cast<TH2F*>(GetESDsData(kPid3));
382       h->Fill(pp,dedx);
383     }
384  
385   }
386
387   // Multiplicity related QA
388   AliESDVZERO     *mltV0 =esd->GetVZEROData();
389   const AliMultiplicity *mltITS=esd->GetMultiplicity();
390   if (mltV0)
391     if (mltITS) {
392        Short_t nv0a=mltV0->GetNbPMV0A();
393        Short_t nv0c=mltV0->GetNbPMV0C();
394        Int_t   nits=mltITS->GetNumberOfTracklets();
395        TH2F *h0=dynamic_cast<TH2F*>(GetESDsData(kMlt0));
396        h0->Fill(nits,nv0a);
397        TH2F *h1=dynamic_cast<TH2F*>(GetESDsData(kMlt1));
398        h1->Fill(nits,nv0c);
399     }
400
401
402   TH1 *tpc=GetESDsData(kTrk2); tpc->Sumw2();
403   TH1 *its=GetESDsData(kTrk3); its->Sumw2();
404   TH1 *tof=GetESDsData(kTrk4); tof->Sumw2();
405   GetESDsData(kTrk5)->Divide(its,tpc,1,1.,"b");
406   GetESDsData(kTrk6)->Divide(tof,tpc,1,1.,"b");
407
408   // V0 related QA
409   Int_t nV0=esd->GetNumberOfV0s();
410   for (Int_t i=0; i<nV0; i++) {
411     Double_t mass;
412     AliESDv0 v0(*esd->GetV0(i));
413
414     v0.ChangeMassHypothesis(kK0Short);
415     mass=v0.GetEffMass();
416     if (v0.GetOnFlyStatus())
417        GetESDsData(kK0on)->Fill(mass);
418     else
419        GetESDsData(kK0off)->Fill(mass);
420
421     v0.ChangeMassHypothesis(kLambda0);
422     mass=v0.GetEffMass();
423     if (v0.GetOnFlyStatus())
424        GetESDsData(kL0on)->Fill(mass);
425     else
426        GetESDsData(kL0off)->Fill(mass);
427
428     v0.ChangeMassHypothesis(kLambda0Bar);
429     mass=v0.GetEffMass();
430     if (v0.GetOnFlyStatus())
431        GetESDsData(kL0on)->Fill(mass);
432     else
433        GetESDsData(kL0off)->Fill(mass);
434   }
435
436 }