]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/hfe/AliAnalysisTaskCheckV0tender.cxx
Fixes to run w/o OCDB infotrmation (in that case the class only provides statistics...
[u/mrichter/AliRoot.git] / PWG3 / hfe / AliAnalysisTaskCheckV0tender.cxx
1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercial purposes is hereby granted   *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  *
14 **************************************************************************/
15 //
16 // Task fir checking the performance of the V0 tender
17 // 
18 // 
19 // Authors
20 //   Matus Kalisky <matus.kalisky@cern.ch>
21 //
22
23 #include <TH1F.h>
24 #include <TList.h>
25
26 #include "AliAnalysisManager.h"
27 #include "AliMCEventHandler.h"
28 #include "AliESDInputHandler.h"
29 #include "AliESDv0.h"
30 #include "AliESDtrack.h"
31 #include "AliMCParticle.h"
32 #include "AliMCEvent.h"
33
34 #include "AliHFEtools.h"
35 #include "AliHFEcollection.h"
36
37 #include "AliAnalysisTaskCheckV0tender.h"
38
39 ClassImp(AliAnalysisTaskCheckV0tender)
40
41 //__________________________________________________________
42 AliAnalysisTaskCheckV0tender::AliAnalysisTaskCheckV0tender():
43   AliAnalysisTaskSE("CheckV0tenderTask")
44   , fOutput(0x0)
45   , fColl(0x0)
46   , fCollMC(0x0)
47   , fEvents(0x0)
48 {
49   //
50   // Default Constructor
51   //
52 }
53 //__________________________________________________________
54 AliAnalysisTaskCheckV0tender::AliAnalysisTaskCheckV0tender(const Char_t *name):
55   AliAnalysisTaskSE(name)
56   , fOutput(0x0)
57   , fColl(0x0)
58   , fCollMC(0x0)
59   , fEvents(0x0)
60 {
61   //
62   // Default Constructor
63   //
64   DefineOutput(1, TH1F::Class());
65   DefineOutput(2, TList::Class());
66
67 }
68 //__________________________________________________________
69 AliAnalysisTaskCheckV0tender::~AliAnalysisTaskCheckV0tender(){
70   //
71   // Destructor
72   //
73
74   if (fOutput) delete fOutput;
75   if (fColl) delete fColl;
76   if (fCollMC) delete fCollMC;
77   if (fEvents) delete fEvents;
78 }
79 //__________________________________________________________
80 void AliAnalysisTaskCheckV0tender::UserCreateOutputObjects(){
81   //
82   // prepare output objects
83   //
84
85   fOutput = new TList();
86   fOutput->SetOwner();
87   // Counter for number of events
88   fEvents = new TH1I("nEvents", "NumberOfEvents", 1, 1, 2);
89
90   fColl = new AliHFEcollection("V0_QA", "tender V0s for data");
91   fCollMC = new AliHFEcollection("V0_MC_QA", "tender V0s for MC");
92
93   // 
94   // Data histos
95   //
96   
97   // total number of tagged V0s
98   fColl->CreateTH1F("h_NumberOf_V0s", "Number of tagged V0s; type; counts", 4, -0.5, 3.5);
99   // pT spectra of the tagged V0s
100   fColl->CreateTH1F("h_Gamma_pt", "p_{T} spectrum of tagged gammas; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
101   fColl->CreateTH1F("h_K0_pt", "p_{T} spectrum of tagged K0s; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
102   fColl->CreateTH1F("h_Lambda_pt", "p_{T} spectrum of tagged Lambdas; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
103   fColl->CreateTH1F("h_ALambda_pt", "p_{T} spectrum of tagged A-Lambdas; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
104   // invariant mass of the V0s
105   fColl->CreateTH1F("h_Gamma_mass", "Inv. Mass of the ", 100, 0., 0.1);
106   fColl->CreateTH1F("h_K0_mass", "Inv. Mass of the ", 100, 0.45, 0.55);
107   fColl->CreateTH1F("h_Lambda_mass", "Inv. Mass of the ", 100, 1.08, 1.14);
108   fColl->CreateTH1F("h_ALambda_mass", "Inv. Mass of the ", 100, 1.08, 1.14);
109
110   // total number of tagged daughter particles (should correlate with number of V0s !
111   fColl->CreateTH1F("h_NumberOfDaughters", "Number of tagged daughters; type; counts", 3, -0.5, 2.5);
112   // pT spectra of tagged daughter particles
113   fColl->CreateTH1F("h_Electron_pt", "p_{T} spectrum of tagged; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
114   fColl->CreateTH1F("h_Pion_pt", "p_{T} spectrum of tagged; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
115   fColl->CreateTH1F("h_Proton_pt", "p_{T} spectrum of tagged; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
116
117   //
118   // MC histos
119   //
120
121   // pT spectra of the tagged V0s
122   fCollMC->CreateTH1F("h_Gamma_pt_S", "MC-S: p_{T} spectrum of tagged gammas; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
123   fCollMC->CreateTH1F("h_K0_pt_S", "MC-S: p_{T} spectrum of tagged K0s; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
124   fCollMC->CreateTH1F("h_Lambda_pt_S", "MC-S: p_{T} spectrum of tagged Lambdas; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
125   fCollMC->CreateTH1F("h_ALambda_pt_S", "MC-S: p_{T} spectrum of tagged A-Lambdas; p_{T} (GeV/c); counts", 20, 0.1, 20, 0); 
126   fCollMC->CreateTH1F("h_Gamma_pt_B", "MC-B: p_{T} spectrum of tagged gammas; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
127   fCollMC->CreateTH1F("h_K0_pt_B", "MC-B: p_{T} spectrum of tagged K0s; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
128   fCollMC->CreateTH1F("h_Lambda_pt_B", "MC-B: p_{T} spectrum of tagged Lambdas; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
129   fCollMC->CreateTH1F("h_ALambda_pt_B", "MC-B: p_{T} spectrum of tagged A-Lambdas; p_{T} (GeV/c); counts", 20, 0.1, 20, 0); 
130   // invariant mass of the V0s
131   fCollMC->CreateTH1F("h_Gamma_mass_S", "MC-S: Inv. Mass of the gamma; m (GeV/c^{2}); counts", 100, 0., 0.1);
132   fCollMC->CreateTH1F("h_K0_mass_S", "MC-S: Inv. Mass of the K0; m (GeV/c^{2}); counts", 100, 0.45, 0.55);
133   fCollMC->CreateTH1F("h_Lambda_mass_S", "MC-S: Inv. Mass of the Lambda; m (GeV/c^{2}); counts", 100, 1.08, 1.14);
134   fCollMC->CreateTH1F("h_ALambda_mass_S", "MC-S: Inv. Mass of the A-Lambda; m (GeV/c^{2}); counts", 100, 1.08, 1.14);
135   fCollMC->CreateTH1F("h_Gamma_mass_B", "MC-B: Inv. Mass of the gamma; m (GeV/c^{2}); counts", 100, 0., 0.1);
136   fCollMC->CreateTH1F("h_K0_mass_B", "MC-B: Inv. Mass of the K0; m (GeV/c^{2}); counts", 100, 0.45, 0.55);
137   fCollMC->CreateTH1F("h_Lambda_mass_B", "MC-B: Inv. Mass of the Lambda; m (GeV/c^{2}); counts", 100, 1.08, 1.14);
138   fCollMC->CreateTH1F("h_ALambda_mass_B", "MC-B: Inv. Mass of the A-Lambda; m (GeV/c^{2}); counts", 100, 1.08, 1.14);
139   // pT spectra of tagged daughter particles
140   fCollMC->CreateTH1F("h_Electron_pt_S", "MC-S: p_{T} spectrum of tagged electrons; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
141   fCollMC->CreateTH1F("h_Pion_pt_S", "MC-S: p_{T} spectrum of tagged pions; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
142   fCollMC->CreateTH1F("h_Proton_pt_S", "MC-S: p_{T} spectrum of tagged protons; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
143   fCollMC->CreateTH1F("h_Electron_pt_B", "MC-B: p_{T} spectrum of tagged electrons; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
144   fCollMC->CreateTH1F("h_Pion_pt_B", "MC-B: p_{T} spectrum of tagged pions; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
145   fCollMC->CreateTH1F("h_Proton_pt_B", "MC-B: p_{T} spectrum of tagged protons; p_{T} (GeV/c); counts", 20, 0.1, 20, 0);
146
147
148   TList *tmp = fColl->GetList();
149   tmp->SetName(fColl->GetName());
150   fOutput->Add(tmp);
151   tmp = 0x0;
152   tmp = fCollMC->GetList();
153   tmp->SetName(fCollMC->GetName());
154   fOutput->Add(tmp);
155   
156   
157 }
158 //__________________________________________________________
159 void AliAnalysisTaskCheckV0tender::UserExec(Option_t *){
160   //
161   // Event Loop
162   // 
163
164   AliMCEventHandler* mcHandler = (dynamic_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler()));
165   //AliESDInputHandler *inh = dynamic_cast<AliESDInputHandler *>(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
166   //AliESDpid *workingPID = NULL;
167   //if(inh && (workingPID = inh->GetESDpid())) workingPID = inh->GetESDpid();
168   //else workingPID = AliHFEtools::GetDefaultPID(mcHandler ? kTRUE : kFALSE);
169      
170   // check the MC data
171   if(fMCEvent && !mcHandler ) return;
172   if(fMCEvent &&  !mcHandler->InitOk() ) return;
173   if(fMCEvent &&  !mcHandler->TreeK() ) return;
174   if(fMCEvent &&  !mcHandler->TreeTR() ) return;
175
176   ProcessV0s();
177   ProcessDaughters();
178
179   if(fMCEvent){
180     ProcessV0sMC();
181     ProcessDaughtersMC();
182   }
183
184   fEvents->Fill(1.1);
185   PostData(1, fEvents);
186   PostData(2, fOutput);
187 }
188 //__________________________________________________________
189 void AliAnalysisTaskCheckV0tender::Terminate(Option_t *){
190   //
191   // Do Post Processing
192   //
193 }
194 //__________________________________________________________
195 void AliAnalysisTaskCheckV0tender::ProcessV0s(){
196   //
197   // loop over the V0s and extract the information about
198   // the V0s tagged by the V0 tender
199   //
200
201   const char * type[4] = {"Gamma", "K0", "Lambda", "ALambda"};
202   char name[256] = "";
203
204   Int_t nV0s = fInputEvent->GetNumberOfV0s();
205   for(Int_t i=0; i<nV0s; ++i){
206     AliESDv0 *esdV0 = (static_cast<AliESDEvent *>(fInputEvent))->GetV0(i);
207     if(!esdV0) continue;
208     if(!esdV0->GetOnFlyStatus()) continue; // Take only V0s from the On-the-fly v0 finder
209     Int_t pid = GetTenderPidV0(esdV0);
210     if(pid < 0) continue;
211     fColl->Fill("h_NumberOf_V0s", pid);
212     sprintf(name, "h_%s_pt", type[pid]);
213     Float_t pT = esdV0->Pt();
214     fColl->Fill(name, pT);
215     Float_t mass = MassV0(esdV0, pid);
216     sprintf(name, "h_%s_mass", type[pid]);
217     fColl->Fill(name, mass);
218   }
219
220 }
221 //__________________________________________________________
222 void AliAnalysisTaskCheckV0tender::ProcessDaughters(){
223   //
224   // produce some check plots for V0 tender tagged single tracks
225   //
226
227   const char * type[3] = {"Electron", "Pion", "Proton"};
228   char name[256] = "";
229
230   Int_t nTracks = fInputEvent->GetNumberOfTracks();
231   for(Int_t i=0; i<nTracks; ++i){
232     AliESDtrack *track = dynamic_cast<AliESDtrack*>(fInputEvent->GetTrack(i));
233     if(!track) continue;
234     Int_t  pid = GetTenderPidDaughter(track);
235     if(pid < 0) continue;
236     fColl->Fill("h_NumberOfDaughters", pid*1.0);
237     Float_t pT = track->Pt();
238     sprintf(name, "h_%s_pt", type[pid]);
239     fColl->Fill(name, pT);
240     
241   }
242 }
243 //__________________________________________________________
244 void AliAnalysisTaskCheckV0tender::ProcessV0sMC(){
245   //
246   // check all V0tender selected V0 on their true identity
247   // 
248
249   const Int_t pid2pdg[4] = {22, 310, 3122, -3122};
250   const char * type[4] = {"Gamma", "K0", "Lambda", "ALambda"};
251   char name[256] = "";
252   Int_t nV0s = fInputEvent->GetNumberOfV0s();
253
254   Int_t nTracks = fInputEvent->GetNumberOfTracks();
255
256   // V0 loop
257   for(Int_t i=0; i<nV0s; ++i){
258     Bool_t id = kFALSE;
259     AliESDv0 *esdV0 = (static_cast<AliESDEvent *>(fInputEvent))->GetV0(i);
260     if(!esdV0) continue;
261     if(!esdV0->GetOnFlyStatus()) continue; // Take only V0s from the On-the-fly v0 finder
262     Int_t pid = GetTenderPidV0(esdV0);
263     if(pid < 0) continue;
264
265     // both ESD daughtr tracks
266     Int_t iN, iP;
267     iN = iP = -1;
268     iP = esdV0->GetPindex();
269     iN = esdV0->GetNindex();    
270     if(iN < 0 || iP < 0) continue;
271     if(iN >= nTracks || iP >= nTracks) continue;    
272     AliESDtrack *dP, *dN;
273     dP = dynamic_cast<AliESDtrack *>(fInputEvent->GetTrack(iP));
274     dN = dynamic_cast<AliESDtrack *>(fInputEvent->GetTrack(iN));  
275     if(!dN || !dP) continue;
276
277     // MC labels of the daughter tracks
278     Int_t lN, lP;
279     lN = dN->GetLabel();
280     lP = dP->GetLabel();
281     if(lN < 0 || lP < 0) continue;
282
283     // MC daughter particles
284     AliMCParticle *mcP, *mcN;
285     mcP = dynamic_cast<AliMCParticle*>(fMCEvent->GetTrack(lP));
286     mcN = dynamic_cast<AliMCParticle*>(fMCEvent->GetTrack(lN));
287     if(!mcP || !mcN) continue;
288
289     // labels of the mother particles
290     Int_t lPm, lNm;
291     lPm = mcP->GetMother();
292     lNm = mcN->GetMother();
293     if(lPm < 0) continue;
294     AliMCParticle *m = dynamic_cast<AliMCParticle*>(fMCEvent->GetTrack(lPm)); 
295     if(!m) continue;
296     Int_t pdg = m->PdgCode();
297     if((lPm == lNm) && (pdg == pid2pdg[pid])) id = kTRUE;
298
299     if(id) sprintf(name, "h_%s_pt_S", type[pid]);
300     else sprintf(name, "h_%s_pt_B", type[pid]);
301     Float_t pT = esdV0->Pt();
302     fCollMC->Fill(name, pT);
303
304     if(id) sprintf(name, "h_%s_mass_S", type[pid]);
305     else sprintf(name, "h_%s_mass_B", type[pid]);
306     Float_t mass = MassV0(esdV0, pid);
307     fCollMC->Fill(name, mass);
308
309    }
310 }
311 //__________________________________________________________
312 void AliAnalysisTaskCheckV0tender::ProcessDaughtersMC(){
313   //
314   // check the identity of the V0tender selected V0 daughters
315   // !!! for positive check only the true identity plays a role here, 
316   // not the true V0 mother identity (e.g. selected electron could come
317   // from primary vertex or pion dalitz deca or true gamma conversion) !!!
318   //
319
320   const Int_t pid2pdg [3] = {11, 211, 2212};
321   const char * type[3] = {"Electron", "Pion", "Proton"};
322   char name[256] = "";
323
324
325   Int_t nTracks = fInputEvent->GetNumberOfTracks();
326   for(Int_t i=0; i<nTracks; ++i){
327     Bool_t id = kFALSE;
328     AliESDtrack *track = dynamic_cast<AliESDtrack*>(fInputEvent->GetTrack(i));
329     if(!track) continue;
330     Int_t  pid = GetTenderPidDaughter(track);
331     if(pid < 0) continue;
332     Float_t pT = track->Pt();
333     Int_t label = track->GetLabel();
334     if(label < 0) continue;
335     AliMCParticle *mcp = 0x0;
336     mcp = dynamic_cast<AliMCParticle*>(fMCEvent->GetTrack(label));
337     if(!mcp) continue;
338     Int_t pdg = TMath::Abs(mcp->PdgCode());
339     if(pdg == pid2pdg[pid]) id = kTRUE;
340     if(id) sprintf(name, "h_%s_pt_S", type[pid]);
341     else sprintf(name, "h_%s_pt_B", type[pid]);
342     fCollMC->Fill(name, pT);
343   }
344 }
345 //__________________________________________________________
346 Int_t AliAnalysisTaskCheckV0tender::GetTenderPidV0(AliESDv0 * const v0){
347   //
348   // retrieve the PID nformation stored in the status flags by the train tender
349   // 
350   Int_t pid = -1;
351   if(!v0){
352     return pid;
353   }
354   Int_t nTimes = 0;
355   if(v0->TestBit(BIT(14))){
356     pid = 0;
357     nTimes++;
358   }
359   if(v0->TestBit(BIT(15))){
360     pid = 1;
361     nTimes++;
362   }
363   if(v0->TestBit(BIT(16))){
364     pid = 2;
365     nTimes++;
366   }
367   if(v0->TestBit(BIT(17))){
368     pid = 3;
369     nTimes++;
370   }
371   if(nTimes > 1){
372     AliWarning("V0 track labeled multiple times by the V0 tender");
373     pid = -1;
374   }    
375
376   //printf(" -D: pid: %i \n", pid);
377
378   return pid;
379 }
380 //__________________________________________________________
381 Int_t AliAnalysisTaskCheckV0tender::GetTenderPidDaughter(AliESDtrack * const track){
382   //
383   // retrieve the PID nformation stored in the status flags by the train tender
384   // 
385
386   Int_t pid = -1;
387   if(!track){
388     return pid;
389   }
390   Int_t nTimes = 0;
391   if(track->TestBit(BIT(14))){
392     pid = 0;
393     nTimes++;
394   }
395   if(track->TestBit(BIT(15))){
396     pid = 1;
397     nTimes++;
398   }
399   if(track->TestBit(BIT(16))){
400     pid = 2;
401     nTimes++;
402   }
403   if(nTimes > 1){
404     AliWarning("V0 track labeled multiple times by the V0 tender");
405     pid = -1;
406   }    
407   return pid;
408 }
409 //__________________________________________________________
410 Float_t AliAnalysisTaskCheckV0tender::MassV0(AliESDv0 * const v0, Int_t id){
411   //
412   // Get the V0 effective mass
413   //
414
415   Float_t mass = -0.1;
416   Bool_t sign = CheckSigns(v0);
417   if(0 == id){
418     mass = v0->GetEffMass(0, 0);
419   }
420   else if(1 == id){
421     mass = v0->GetEffMass(2, 2);
422   }
423   else if(2 == id){
424     mass = (sign) ? v0->GetEffMass(4, 2) : v0->GetEffMass(2, 4);
425   }
426   else if(3 == id){
427     mass = (sign) ? v0->GetEffMass(2, 4) : v0->GetEffMass(4, 2);
428   }
429   else{
430     AliWarning(Form("Unrecognized V0 id: %i", id));
431   }
432
433   return mass;
434
435 }
436 //__________________________________________________________
437 Bool_t AliAnalysisTaskCheckV0tender::CheckSigns(AliESDv0 * const v0){
438   //
439   // check wheter the sign was correctly applied to 
440   // V0 daughter tracks
441   // This function should become obsolete once the V0 finder will be updated (fixed)
442   //
443   
444   Bool_t correct = kFALSE;
445
446   Int_t pIndex = 0, nIndex = 0;
447   pIndex = v0->GetPindex();
448   nIndex = v0->GetNindex();
449   
450   AliESDtrack* d[2];
451   d[0] = dynamic_cast<AliESDtrack*>(fInputEvent->GetTrack(pIndex));
452   d[1] = dynamic_cast<AliESDtrack*>(fInputEvent->GetTrack(nIndex));
453
454   Int_t sign[2];
455   sign[0] = (int)d[0]->GetSign();
456   sign[1] = (int)d[1]->GetSign();
457   
458   if(-1 == sign[0] && 1 == sign[1]){
459     correct = kFALSE;
460     //v0->SetIndex(0, pIndex);  // set the index of the negative v0 track
461     //v0->SetIndex(1, nIndex);  // set the index of the positive v0 track
462   }
463   else{
464     correct = kTRUE;
465   }
466   
467   //pIndex = v0->GetPindex();
468   //nIndex = v0->GetNindex();
469   //printf("-D2: P: %i, N: %i\n", pIndex, nIndex);
470
471   return correct;
472 }