AliPhysics  b095172 (b095172)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliAnalysisTaskPtEfficiencyJets.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * Copyright(c) 1998-2015, 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 #include <TClonesArray.h>
16 #include <TList.h>
17 #include <TMath.h>
18 #include <TNtuple.h>
19 
21 #include "AliAnalysisUtils.h"
22 #include "AliAODMCParticle.h"
23 #include "AliEmcalJet.h"
24 #include "AliGenPythiaEventHeader.h"
25 #include "AliJetContainer.h"
26 #include "AliMCEvent.h"
27 #include "AliMCParticle.h"
28 #include "AliParticleContainer.h"
29 #include "AliVEvent.h"
30 #include "AliVParticle.h"
31 #include "AliVTrack.h"
32 #include "AliVVertex.h"
33 #include "AliEmcalTrackSelection.h"
34 
36 
38 
39 namespace EMCalTriggerPtAnalysis {
40 
41 AliAnalysisTaskPtEfficiencyJets::AliAnalysisTaskPtEfficiencyJets() :
43  fAnalysisUtils(NULL),
44  fMCJetContainer(""),
45  fTrackCuts(NULL),
46  fTrackNtuple(NULL)
47 {
48 }
49 
50 AliAnalysisTaskPtEfficiencyJets::AliAnalysisTaskPtEfficiencyJets(const char *name) :
51  AliAnalysisTaskEmcalJet(name, kTRUE),
52  fAnalysisUtils(NULL),
53  fMCJetContainer(""),
54  fTrackCuts(NULL),
55  fTrackNtuple(NULL)
56 {
57  fAnalysisUtils = new AliAnalysisUtils();
58  SetMakeGeneralHistograms(kTRUE);
59 }
60 
61 AliAnalysisTaskPtEfficiencyJets::~AliAnalysisTaskPtEfficiencyJets() {
62  if(fAnalysisUtils) delete fAnalysisUtils;
63  if(fTrackCuts) delete fTrackCuts;
64  if(fTrackNtuple) delete fTrackNtuple;
65 }
66 
67 void AliAnalysisTaskPtEfficiencyJets::UserCreateOutputObjects() {
68  OpenFile(1);
70  fTrackNtuple = new TNtuple("tracktuple", "track tuple", "partpt:trackpt:jetpt:pthard:parteta:partphi:tracketa:trackphi:vertexz");
71  fOutput->Add(fTrackNtuple);
72  PostData(1, fOutput);
73 }
74 
75 Bool_t AliAnalysisTaskPtEfficiencyJets::Run() {
76  // Apply event selection first
77  if(!fAnalysisUtils->IsVertexSelected2013pA(fInputEvent)) return kFALSE;
78  const AliVVertex *primvertex = fInputEvent->GetPrimaryVertex();
79 
80  Float_t pthard = 0;
81  AliGenPythiaEventHeader *evheader = dynamic_cast<AliGenPythiaEventHeader *>(fMCEvent->GenEventHeader());
82  if(evheader) pthard = evheader->GetPtHard();
83  // get the Monte-Carlo jet container
84  AliJetContainer *jcmc = dynamic_cast<AliJetContainer *>(fJetCollArray.FindObject(fMCJetContainer.Data()));
85  if(!jcmc) return kFALSE;
86  TClonesArray *particleArray = jcmc->GetParticleContainer()->GetArray();
87 
88  for(TIter partIter = TIter(particleArray).Begin(); partIter != TIter::End(); ++partIter){
89  AliVParticle *part = static_cast<AliVParticle *>(*partIter);
90  if(!SelectTrueParticle(part)) continue;
91  AliVTrack *reconstructed = FindAssociatedTrack(part);
92  AliEmcalJet *jet = FindAssociatedJet(part, jcmc);
93 
94  // Fill ntuple with
95  // particle pt
96  // track pt
97  // jet pt
98  // particle eta
99  // particle phi
100  // track eta
101  // track phi
102  // vertex z
103  Float_t data[] = {(Float_t)TMath::Abs(part->Pt()), reconstructed ? static_cast<Float_t>(TMath::Abs(reconstructed->Pt())) : 0,
104  jet ? static_cast<Float_t>(TMath::Abs(jet->Pt())) : 0, pthard, (Float_t)part->Eta(), (Float_t)part->Phi(), static_cast<Float_t>(reconstructed ? reconstructed->Eta() : -1000.),
105  static_cast<Float_t>(reconstructed ? reconstructed->Phi() : -1000.), static_cast<Float_t>(primvertex ? primvertex->GetZ() : -1000.)
106  };
107  fTrackNtuple->Fill(data);
108  }
109  PostData(1, fOutput);
110  return kTRUE;
111 }
112 
113 AliVTrack* AliAnalysisTaskPtEfficiencyJets::FindAssociatedTrack(AliVParticle* trueParticle) {
114  AliVTrack *result = NULL, *tmp = NULL;
115  for(int itrk = 0; itrk < fInputEvent->GetNumberOfTracks(); itrk++){
116  AliVTrack *trk = static_cast<AliVTrack *>(fInputEvent->GetTrack(itrk));
117  if(TMath::Abs(trk->GetLabel()) != TMath::Abs(trueParticle->GetLabel())) continue;
118  if(!fTrackCuts->IsTrackAccepted(trk)) continue;
119  result = trk;
120  break;
121  }
122  return result;
123 }
124 
125 AliEmcalJet* AliAnalysisTaskPtEfficiencyJets::FindAssociatedJet(AliVParticle* trueParticle, AliJetContainer* jets) {
126  AliEmcalJet *result = NULL;
127  for(int ijet = 0; ijet < jets->GetNJets(); ijet++){
128  AliEmcalJet *nextjet = jets->GetJet(ijet);
129  if(nextjet->ContainsTrack(trueParticle, jets->GetParticleContainer()->GetArray())){
130  result = nextjet;
131  break;
132  }
133  }
134  return result;
135 }
136 
137 bool AliAnalysisTaskPtEfficiencyJets::SelectTrueParticle(AliVParticle* part) {
138  if(part->IsA() == AliAODMCParticle::Class()){
139  AliAODMCParticle *aodpart = static_cast<AliAODMCParticle *>(part);
140  if(!aodpart->IsPhysicalPrimary()) return kFALSE;
141  if(!aodpart->Charge()) return kFALSE;
142  if(TMath::Abs(aodpart->Eta()) > 0.8) return kFALSE;
143  if(TMath::Abs(aodpart->Pt()) < 5.) return kFALSE;
144  return kTRUE;
145  } else {
146  AliMCParticle *esdpart = static_cast<AliMCParticle *>(part);
147  if(!fMCEvent->IsPhysicalPrimary(esdpart->GetLabel())) return kFALSE;
148  if(!esdpart->Charge()) return kFALSE;
149  if(TMath::Abs(esdpart->Eta()) > 0.8) return kFALSE;
150  if(TMath::Abs(esdpart->Pt()) < 5.) return kFALSE;
151  return kTRUE;
152  }
153 }
154 
155 } /* namespace EMCalTriggerPtAnalysis */
AliParticleContainer * GetParticleContainer() const
float Float_t
Definition: External.C:68
Int_t GetNJets() const
Int_t ContainsTrack(AliVParticle *track, TClonesArray *tracks) const
Double_t Pt() const
Definition: AliEmcalJet.h:93
Base task in the EMCAL jet framework.
Represent a jet reconstructed using the EMCal jet framework.
Definition: AliEmcalJet.h:44
bool Bool_t
Definition: External.C:53
ClassImp(EMCalTriggerPtAnalysis::AliAnalysisTaskPtEfficiencyJets) namespace EMCalTriggerPtAnalysis
Container for jet within the EMCAL jet framework.
TList * OpenFile(const char *fname)
Definition: DrawAnaELoss.C:65
AliEmcalJet * GetJet(Int_t i) const