]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/Base/AliAnalysisTaskCFTree.cxx
fc6e4ec9935d29b1453fff69aa594288603b3f47
[u/mrichter/AliRoot.git] / PWGCF / Correlations / Base / AliAnalysisTaskCFTree.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 // Analysis task to produce trees of lightweight events
17 // evgeny.kryshen@cern.ch
18
19 // aliroot
20 #include "AliAnalysisTaskCFTree.h"
21 #include "AliCFParticle.h"
22 #include "AliAnalysisManager.h"
23 #include "AliInputEventHandler.h"
24 #include "AliESDEvent.h"
25 #include "AliAODEvent.h"
26 #include "AliMCEvent.h"
27 #include "AliHeader.h"
28 #include "AliESDHeader.h"
29 #include "AliAODHeader.h"
30 #include "AliAODMCHeader.h"
31 #include "AliGenEventHeader.h"
32 #include "AliGenCocktailEventHeader.h"
33 #include "AliAODTrack.h"
34 #include "AliESDtrack.h"
35 #include "AliMCParticle.h"
36 #include "AliAODMCParticle.h"
37 #include "AliExternalTrackParam.h"
38 #include "AliCentrality.h"
39 #include "AliAnalysisFilter.h"
40 #include "AliPIDResponse.h"
41 #include "AliTPCPIDResponse.h"
42 #include "AliAODTracklets.h"
43 // root
44 #include "TMath.h"
45 #include "TFile.h"
46 #include "TList.h"
47 #include "TH1I.h"
48 #include "TChain.h"
49 #include "TTree.h"
50 #include "TClonesArray.h"
51 ClassImp(AliAnalysisTaskCFTree)
52
53 //-----------------------------------------------------------------------------
54 AliAnalysisTaskCFTree::AliAnalysisTaskCFTree(const char* name) :
55   AliAnalysisTask(name,""),
56   fDebug(0),
57   fMode(0),
58   fIsAOD(1),
59   fInputHandler(0x0),
60   fMcHandler(0x0),
61   fTrackFilter(0x0),
62   fHybridConstrainedMask(0),
63   fTPConlyConstrainedMask(0),
64   fPIDResponse(0x0),
65   fListOfHistos(0x0),
66   fEventStatistics(0x0),
67   fTree(0x0),
68   fParticles(0x0),
69   fTracklets(0x0),
70   fMuons(0x0),
71   fField(0),
72   fCentrality(0),
73   fZvtx(0),
74   fRunNumber(0),
75   fPeriod(0),
76   fOrbit(),
77   fBc(),
78   fSelectMask(0),
79   fSelectBit(AliVEvent::kMB),
80   fZVertexCut(10.),
81   fnTracksVertex(1),
82   fCentralityMethod("V0M"),
83   fTrackFilterBit(0xffffffff),
84   fTrackEtaCut(1.0),
85   fPtMin(0.15),
86   fSharedClusterCut(0.4),
87   fCrossedRowsCut(100),
88   fFoundFractionCut(0.8),
89   fStoreTracks(0),
90   fStoreTracklets(0),
91   fStoreMuons(0)
92 {
93   Info("AliAnalysisTaskCFTree","Calling Constructor");
94   DefineInput(0,TChain::Class());
95   DefineOutput(0,TList::Class());
96   DefineOutput(1,TTree::Class());
97 }
98 //-----------------------------------------------------------------------------
99
100
101 //-----------------------------------------------------------------------------
102 void AliAnalysisTaskCFTree::ConnectInputData(Option_t* /*option*/){
103   fInputHandler = dynamic_cast<AliInputEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
104   fMcHandler    = dynamic_cast<AliInputEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
105   fInputHandler->SetNeedField();
106 }
107 //-----------------------------------------------------------------------------
108
109
110 //-----------------------------------------------------------------------------
111 void AliAnalysisTaskCFTree::CreateOutputObjects(){
112   fListOfHistos = new TList();
113   fListOfHistos->SetOwner();
114   fEventStatistics = new TH1I("fEventStatistics","",10,0,10);
115   fEventStatistics->SetBit(TH1::kCanRebin);
116
117   fListOfHistos->Add(fEventStatistics);
118
119   if (fStoreTracks)    fParticles = new TClonesArray("AliCFParticle",2000);
120   if (fStoreTracklets) fTracklets = new TClonesArray("AliCFParticle",100);
121   if (fStoreMuons)     fMuons     = new TClonesArray("AliCFParticle",100);
122   // create file-resident tree
123   TDirectory *owd = gDirectory;
124   OpenFile(1);
125   fTree = new TTree("events","events");
126   owd->cd();
127   fTree->Branch("cent",&fCentrality);
128   fTree->Branch("zvtx",&fZvtx);
129   fTree->Branch("field",&fField);
130   fTree->Branch("run",&fRunNumber);
131   fTree->Branch("period",&fPeriod);
132   fTree->Branch("orbit",&fOrbit);
133   fTree->Branch("bc",&fBc);
134   fTree->Branch("mask",&fSelectMask);
135   if (fStoreTracks)    fTree->Branch("particles",&fParticles);
136   if (fStoreTracklets) fTree->Branch("tracklets",&fTracklets);
137   if (fStoreMuons)     fTree->Branch("muons",&fMuons);
138   PostData(0,fListOfHistos);
139   PostData(1,fTree);
140 }
141 //-----------------------------------------------------------------------------
142
143
144 //-----------------------------------------------------------------------------
145 void AliAnalysisTaskCFTree::Exec(Option_t *){
146   if (fParticles) fParticles->Clear();
147   if (fTracklets) fTracklets->Clear();
148   if (fMuons)     fMuons->Clear();
149   
150   fEventStatistics->Fill("before cuts",1);
151   
152   AliVEvent* event = fInputHandler->GetEvent();
153   if (!event) return;
154   fEventStatistics->Fill("after event check",1);
155   
156   fPIDResponse = fInputHandler->GetPIDResponse();
157   
158   fSelectMask = fInputHandler->IsEventSelected();
159   if (!(fSelectMask & fSelectBit)) return;
160   fEventStatistics->Fill("after trigger selection",1);
161   fRunNumber  = event->GetRunNumber();
162   fPeriod     = event->GetPeriodNumber();
163   fOrbit      = event->GetOrbitNumber();
164   fBc         = event->GetBunchCrossNumber();
165   fField      = event->GetMagneticField();
166   fCentrality = event->GetCentrality()->GetCentralityPercentile(fCentralityMethod);
167   
168   if (fMode==0){ // Data analysis
169     const AliVVertex* vertex  = event->GetPrimaryVertex();
170     if (!vertex) return;
171     fZvtx  = vertex->GetZ();
172     TString name(vertex->GetName());
173     if (name.CompareTo("PrimaryVertex") && name.CompareTo("SPDVertex")) return; // Reject TPC only vertex
174     fEventStatistics->Fill("after rejection of TPC only vertex",1);
175     if (TMath::Abs(fZvtx) >= fZVertexCut || vertex->GetNContributors()<fnTracksVertex)  return;
176     fEventStatistics->Fill("after vertex selection",1);
177     
178     if (fParticles) {
179       for (Int_t ipart=0;ipart<event->GetNumberOfTracks();ipart++){
180         AliVTrack* track = (AliVTrack*) event->GetTrack(ipart);
181         if (!track) continue;
182         UInt_t mask = GetFilterMap(track);
183         if (!(mask & fTrackFilterBit)) continue;
184
185         if (track->InheritsFrom("AliAODTrack")) AddTrack(track,mask,0);
186         else if (track->InheritsFrom("AliESDtrack")) {
187           if (mask)                           AddTrack(track,mask,1);
188           if (mask & fHybridConstrainedMask)  AddTrack(track,mask,2);
189           if (mask & fTPConlyConstrainedMask) AddTrack(track,mask,3);
190         }
191       }
192     }
193     if (fIsAOD){
194       AliAODEvent* aod = (AliAODEvent*) event;
195       
196       if (fStoreTracklets){
197         AliAODTracklets* tracklets = aod->GetTracklets();
198         Int_t nTracklets = tracklets->GetNumberOfTracklets();
199         for (Int_t i=0;i<nTracklets;i++){
200           Float_t phi   = tracklets->GetPhi(i);
201           Float_t theta = tracklets->GetTheta(i);
202           Float_t dphi  = tracklets->GetDeltaPhi(i);
203           new ((*fTracklets)[fTracklets->GetEntriesFast()]) AliCFParticle(dphi,-TMath::Log(TMath::Tan(theta/2)),phi,0,0);
204         }
205       }
206       
207       if (fStoreMuons){
208         for (Int_t iTrack = 0; iTrack < aod->GetNTracks(); iTrack++) {
209           AliAODTrack* track = aod->GetTrack(iTrack);
210           if (!track->IsMuonTrack()) continue;
211           Float_t pt     = track->Pt();
212           Float_t eta    = track->Eta();
213           Float_t phi    = track->Phi();
214           Short_t charge = track->Charge();
215           Float_t dca    = track->DCA();
216           Float_t chi2   = track->Chi2perNDF();
217           Float_t rabs   = track->GetRAtAbsorberEnd();
218           Int_t   mask   = track->GetMatchTrigger();
219           if (rabs < 17.6 || rabs > 89.5) continue;
220           if (eta < -4 || eta > -2.5) continue;
221           AliCFParticle* part = new ((*fMuons)[fMuons->GetEntriesFast()]) AliCFParticle(pt,eta,phi,charge,mask,3);
222           part->SetAt(dca,0);
223           part->SetAt(chi2,1);
224           part->SetAt(rabs,2);
225         }
226       }
227     }
228   }
229   
230   else { // MC analysis
231     AliMCEvent* mcEvent = fMcHandler ? fMcHandler->MCEvent() : 0;
232     TClonesArray* mcTracks = dynamic_cast<TClonesArray*>(event->FindListObject(AliAODMCParticle::StdBranchName()));
233     if (!mcEvent && !mcTracks) { printf("No mc object found\n"); return; }
234     fEventStatistics->Fill("after mc objects check",1);
235
236     Int_t nPrimGen = 0;
237     Int_t nProduced = 0;
238     if (mcEvent) { // ESD
239       AliHeader* header = (AliHeader*) mcEvent->Header();
240       AliGenCocktailEventHeader* cocktailHeader = dynamic_cast<AliGenCocktailEventHeader*> (header->GenEventHeader());
241       AliGenEventHeader* mcHeader = dynamic_cast<AliGenEventHeader*> (cocktailHeader ? cocktailHeader->GetHeaders()->First() : header->GenEventHeader());
242       if (!mcHeader) { printf("mc header not found\n"); return; }
243       nProduced = mcEvent->GetNumberOfTracks();
244       nPrimGen = mcHeader->NProduced();
245       fZvtx = mcEvent->GetPrimaryVertex()->GetZ();
246     } else { // AOD
247       AliAODMCHeader* mcHeader = (AliAODMCHeader*) event->FindListObject(AliAODMCHeader::StdBranchName());
248       if (!mcHeader) { printf("AliAODMCHeader not found\n"); return; }
249       if (mcHeader->GetCocktailHeaders()) {
250         nProduced = mcTracks->GetEntriesFast();
251         AliGenEventHeader* header0 =  mcHeader->GetCocktailHeader(0);
252         if (!header0) { printf("first header expected but not found\n"); return; }
253         nPrimGen = header0->NProduced();
254       } else  nPrimGen = nProduced;
255       fZvtx = mcHeader->GetVtxZ();
256     }
257     fEventStatistics->Fill("after mc header check",1);
258
259     if (TMath::Abs(fZvtx)>fZVertexCut) return;
260     fEventStatistics->Fill("after MC vertex cut",1);
261
262     const AliVVertex* vertex = event->GetPrimaryVertex();
263     if (!vertex) return;
264     fEventStatistics->Fill("after check for vertex object",1);
265     TString name(vertex->GetName());
266     if (name.CompareTo("PrimaryVertex") && name.CompareTo("SPDVertex")) fZvtx=1000; // Reject TPC only vertex
267     fEventStatistics->Fill("after rejection of TPC only vertex",1);
268     if (vertex->GetNContributors()<fnTracksVertex)  fZvtx=-1000;
269     fEventStatistics->Fill("after check for vertex contributors",1);
270     
271     if (fParticles) {
272       UInt_t* masks = new UInt_t[nProduced];
273       for (Int_t i=0;i<nProduced;i++) masks[i]=0;
274       
275       // Loop over reconstructed tracks to set masks
276       for (Int_t ipart=0;ipart<event->GetNumberOfTracks();ipart++){
277         AliVTrack* part = (AliVTrack*) event->GetTrack(ipart);
278         if (!part) continue;
279         Int_t label = TMath::Abs(part->GetLabel());
280         if (label>=nProduced) continue;
281         masks[label] |= GetFilterMap(part);
282       }
283       
284       // Loop over mc tracks to be stored
285       for (Int_t ipart=0;ipart<nProduced;ipart++){
286         AliVParticle* part = 0;
287         Bool_t isPrimary = 0;
288         if (mcTracks) { // AOD analysis
289           AliAODMCParticle* particle = (AliAODMCParticle*) mcTracks->At(ipart);
290           if (!particle) continue;
291           // skip injected signals
292           AliAODMCParticle* mother = particle;
293           while (!mother->IsPrimary()) mother = (AliAODMCParticle*) mcTracks->At(mother->GetMother());
294           if (mother->GetLabel()>=nPrimGen) continue;
295           part = particle;
296           // check for primary
297           isPrimary = particle->IsPhysicalPrimary();
298         } else { // ESD analysis
299           AliMCParticle* particle = (AliMCParticle*) mcEvent->GetTrack(ipart);
300           if (!particle) continue;
301           // skip injected signals
302           AliMCParticle* mother = particle;
303           while (mother->GetMother()>=0) mother = (AliMCParticle*) mcEvent->GetTrack(mother->GetMother());
304           if (mother->GetLabel()>=nPrimGen) continue;
305           part = particle;
306           // check for primary
307           isPrimary = mcEvent->IsPhysicalPrimary(ipart);
308         }
309         
310         // store only primaries and all reconstructed (non-zero mask)
311         Int_t mask = masks[ipart];
312         if (isPrimary) mask |= (1 << 30);
313         if (!mask) continue; 
314         AddTrack(part,mask);
315       }
316       delete[] masks;
317     }
318     
319     if (fIsAOD){
320       AliAODEvent* aod = (AliAODEvent*) event;
321     
322       if (fStoreTracklets) {
323         AliAODTracklets* tracklets = aod->GetTracklets();
324         Int_t nTracklets = tracklets->GetNumberOfTracklets();
325         for (Int_t i=0;i<nTracklets;i++){
326           Float_t phi   = tracklets->GetPhi(i);
327           Float_t theta = tracklets->GetTheta(i);
328           Float_t dphi  = tracklets->GetDeltaPhi(i);
329           AliCFParticle* tracklet = new ((*fTracklets)[fTracklets->GetEntriesFast()]) AliCFParticle(dphi,-TMath::Log(TMath::Tan(theta/2)),phi,0,0,4);
330           Int_t label1 = tracklets->GetLabel(i,0);
331           Int_t label2 = tracklets->GetLabel(i,1);
332           if (label1!=label2) continue;
333           AliAODMCParticle* particle = (AliAODMCParticle*) mcTracks->At(label1);
334           if (!particle) continue;
335           Short_t charge = particle->Charge();
336           Float_t ptMC   = particle->Pt();
337           Float_t etaMC  = particle->Eta();
338           Float_t phiMC  = particle->Phi();
339           Float_t pdg    = particle->PdgCode();
340           tracklet->SetCharge(charge);
341           tracklet->SetAt(ptMC,0);
342           tracklet->SetAt(etaMC,1);
343           tracklet->SetAt(phiMC,2);
344           tracklet->SetAt(pdg,3);
345         }
346       }
347
348       if (fStoreMuons){
349         for (Int_t iTrack = 0; iTrack < aod->GetNTracks(); iTrack++) {
350           AliAODTrack* track = aod->GetTrack(iTrack);
351           if (!track->IsMuonTrack()) continue;
352           Float_t pt     = track->Pt();
353           Float_t eta    = track->Eta();
354           Float_t phi    = track->Phi();
355           Short_t charge = track->Charge();
356           Float_t dca    = track->DCA();
357           Float_t chi2   = track->Chi2perNDF();
358           Float_t rabs   = track->GetRAtAbsorberEnd();
359           Int_t   mask   = track->GetMatchTrigger();
360           if (rabs < 17.6 || rabs > 89.5) continue;
361           if (eta < -4 || eta > -2.5) continue;
362           AliCFParticle* part = new ((*fMuons)[fMuons->GetEntriesFast()]) AliCFParticle(pt,eta,phi,charge,mask,3);
363           part->SetAt(dca,0);
364           part->SetAt(chi2,1);
365           part->SetAt(rabs,2);
366         }
367       }
368     }
369   }
370   fTree->Fill();
371   PostData(0,fListOfHistos);
372   PostData(1,fTree);
373 }
374 //-----------------------------------------------------------------------------
375
376
377 //-----------------------------------------------------------------------------
378 UInt_t AliAnalysisTaskCFTree::GetFilterMap(AliVParticle* particle){
379   UInt_t mask = 0;
380   if (particle->InheritsFrom("AliAODTrack")) {
381     AliAODTrack* part = (AliAODTrack*) particle;
382     mask = part->GetFilterMap();
383     Double_t nCrossedRaws      = part->GetTPCNCrossedRows();
384     Double_t nFindableClusters = part->GetTPCNclsF();
385     Double_t nSharedClusters   = part->GetTPCnclsS();
386     Double_t nClusters         = part->GetTPCncls();
387     Bool_t itsRefit            = part->GetStatus() & AliVTrack::kITSrefit;
388     if (nCrossedRaws/nFindableClusters > fFoundFractionCut) mask |= (1 << 26);
389     if (nCrossedRaws>fCrossedRowsCut)                       mask |= (1 << 27);
390     if (itsRefit)                                           mask |= (1 << 28);
391     if (nSharedClusters/nClusters<=fSharedClusterCut)       mask |= (1 << 29);
392     if (part->GetLabel()<0)                                 mask |= (1 << 31);
393   } else if (particle->InheritsFrom("AliESDtrack")){
394     AliESDtrack* part = (AliESDtrack*) particle;
395     if (!fTrackFilter) AliFatal("Track filter undefined");
396     mask |= fTrackFilter->IsSelected(part);
397   }
398   
399   return mask;
400 }
401 //-----------------------------------------------------------------------------
402
403
404 //-----------------------------------------------------------------------------
405 AliCFParticle* AliAnalysisTaskCFTree::AddTrack(AliVParticle* part, UInt_t mask, UInt_t flag){
406
407   // skip neutral mc particles
408   Char_t charge = part->Charge();
409   if (charge==0) return NULL;
410   
411   // set pt,eta,phi
412   Float_t pt=0,eta=0,phi=0;
413   if (flag==0 || flag==1){ // AOD, MC or Global ESD tracks
414     pt  = part->Pt();
415     eta = part->Eta();
416     phi = part->Phi();
417     if  (flag==1) mask &= (~fHybridConstrainedMask) & (~fTPConlyConstrainedMask);
418   } 
419   else if (flag==2) { // Hybrid constrained tracks (ESD)
420     AliESDtrack* track = (AliESDtrack*) part;
421     const AliExternalTrackParam* param = track->GetConstrainedParam();
422     pt  = param->Pt();
423     eta = param->Eta();
424     phi = param->Phi();
425     mask &= fHybridConstrainedMask;
426   } 
427   else if (flag==3) { // TPC only constrained tracks (ESD)
428     AliESDtrack* track = (AliESDtrack*) part;
429     AliESDtrack tpcTrack;
430     if (!track->FillTPCOnlyTrack(tpcTrack)) return NULL;
431     AliExternalTrackParam param;
432     const AliESDVertex* vtxSPD = ((AliESDEvent*) fInputHandler->GetEvent())->GetPrimaryVertexSPD();
433     if (!tpcTrack.RelateToVertexTPC(vtxSPD,fField,1.e30,&param)) return NULL;
434     pt  = param.Pt();
435     eta = param.Eta();
436     phi = param.Phi();
437     mask &= fTPConlyConstrainedMask;
438   }
439
440   // kinematic cuts
441   if (pt < fPtMin || TMath::Abs(eta) > fTrackEtaCut) return NULL;
442
443   
444   AliCFParticle* cfpart = new ((*fParticles)[fParticles->GetEntriesFast()]) AliCFParticle(pt,eta,phi,charge,mask,3);
445
446   AliVTrack* track = dynamic_cast<AliVTrack*> (part);
447   if (!track) return cfpart;
448   Float_t ncl  = track->GetTPCsignalN();
449   Float_t dedx = track->GetTPCsignalTunedOnData(); if (dedx<=0) dedx = track->GetTPCsignal();
450   Float_t beta = -99;
451   if (track->GetStatus()&AliESDtrack::kTOFpid){
452     Double_t tof[5];
453     track->GetIntegratedTimes(tof);
454     beta = tof[0]/track->GetTOFsignal();
455   }
456   cfpart->SetAt(ncl,0);
457   cfpart->SetAt(dedx,1);
458   cfpart->SetAt(beta,2);
459   return cfpart;
460 }