AliPhysics  eb0e5d9 (eb0e5d9)
AliEmcalJet.h
Go to the documentation of this file.
1 #ifndef ALIEMCALJET_H
2 #define ALIEMCALJET_H
3 /* Copyright(c) 1998-2016, ALICE Experiment at CERN, All rights reserved. *
4  * See cxx source for full Copyright notice */
5 
6 #include <vector>
7 #include <algorithm>
8 #include <utility>
9 
10 #include <iosfwd>
11 #include <TArrayI.h>
12 #include <TMath.h>
13 #include <TClonesArray.h>
14 #include <TVector2.h>
15 #include <TLorentzVector.h>
16 #include <TString.h>
17 
18 #include <AliVParticle.h>
19 #include <AliVCluster.h>
20 #include <AliVEvent.h>
21 
25 
51 class AliEmcalJet : public AliVParticle
52 {
53  public:
54 
67  kTPC = 1<<0,
68  kTPCfid = 1<<1,
69  kEMCAL = 1<<2,
70  kEMCALfid = 1<<3,
71  kDCAL = 1<<4,
72  kDCALfid = 1<<5,
73  kDCALonly = 1<<6,
74  kDCALonlyfid = 1<<7,
75  kPHOS = 1<<8,
76  kPHOSfid = 1<<9,
77  kUser = 1<<10
78  };
79 
84  enum EFlavourTag {
85  kDStar = 1<<0,
86  kD0 = 1<<1,
87  kSig1 = 1<<2,
88  kSig2 = 1<<3,
89  kBckgrd1 = 1<<4,
90  kBckgrd2 = 1<<5,
91  kBckgrd3 = 1<<6
92  };
93 
94  AliEmcalJet();
97  AliEmcalJet(const AliEmcalJet &jet);
98  AliEmcalJet& operator=(const AliEmcalJet &jet);
99  virtual ~AliEmcalJet();
100  friend std::ostream &operator<<(std::ostream &in, const AliEmcalJet &myjet);
101  Int_t Compare(const TObject* obj) const;
102  std::ostream &Print(std::ostream &in) const;
103  TString toString() const;
104 
105  // Implementation of AliVParticle interface
106  Double_t Px() const { return fPt*TMath::Cos(fPhi) ; }
107  Double_t Py() const { return fPt*TMath::Sin(fPhi) ; }
108  Double_t Pz() const { return fPt*TMath::SinH(fEta); }
109  Double_t Pt() const { return fPt ; }
110  Double_t P() const { return fPt*TMath::CosH(fEta); }
111  Bool_t PxPyPz(Double_t p[3]) const { p[0]=Px();p[1]=Py();p[2]=Pz(); return kTRUE; }
112  Double_t Xv() const { return 0.; }
113  Double_t Yv() const { return 0.; }
114  Double_t Zv() const { return 0.; }
115  Bool_t XvYvZv(Double_t x[3]) const { x[0]=0;x[1]=0;x[2]=0 ; return kTRUE; }
116  Double_t OneOverPt() const { return 1./fPt ; }
117  Double_t Phi() const { return fPhi ; }
118  Double_t Theta() const { return 2*TMath::ATan(TMath::Exp(-fEta)) ; }
119  Double_t E() const { Double_t p=P(); return TMath::Sqrt(fM*fM+p*p); }
120  Double_t M() const { return fM ; }
121  Double_t Eta() const { return fEta ; }
122  Double_t Y() const { Double_t e = E(); Double_t pz = Pz(); return 0.5*TMath::Log((e+pz)/(e-pz)); }
123  Short_t Charge() const { return 0 ; }
124  Int_t GetLabel() const { return fLabel ; }
125  Int_t PdgCode() const { return 0; }
126  const Double_t *PID() const { return 0; }
127 
128  // Other kinematic and jet properties
129  Double_t Phi_0_2pi() const { return TVector2::Phi_0_2pi(fPhi); }
130  Double_t Area() const { return fArea ; }
131  Double_t AreaPt() const { return fArea ; }
132  Double_t AreaEta() const { return fAreaEta ; }
133  Double_t AreaPhi() const { return fAreaPhi ; }
134  Double_t AreaE() const { return fAreaE ; }
135  Double_t AreaEmc() const { return fAreaEmc ; }
136  Bool_t AxisInEmcal() const { return fAxisInEmcal ; }
137  Int_t ClusterAt(Int_t idx) const { return fClusterIDs.At(idx) ; }
138  UShort_t GetNumberOfClusters() const { return fClusterIDs.GetSize() ; }
139  UShort_t GetNumberOfTracks() const { return fTrackIDs.GetSize() ; }
141  Double_t FracEmcalArea() const { return fAreaEmc/fArea ; }
142  Bool_t IsInsideEmcal() const { return (fAreaEmc/fArea>0.999) ; }
143  Bool_t IsInEmcal() const { return (Bool_t)(fAreaEmc > 0) ; }
144  Bool_t IsMC() const { return (Bool_t)(MCPt() > 0) ; }
145  Bool_t IsSortable() const { return kTRUE ; }
146  Double_t MaxNeutralPt() const { return fMaxNPt ; }
147  Double_t MaxChargedPt() const { return fMaxCPt ; }
148  Double_t NEF() const { return fNEF ; }
149  UShort_t Nn() const { return fNn ; }
150  UShort_t Nch() const { return fNch ; }
151  UShort_t N() const { return Nch()+Nn() ; }
152  Int_t NEmc() const { return fNEmc ; }
153  Double_t MCPt() const { return fMCPt ; }
154  Double_t MaxClusterPt() const { return MaxNeutralPt() ; }
155  Double_t MaxTrackPt() const { return MaxChargedPt() ; }
156  Double_t MaxPartPt() const { return fMaxCPt < fMaxNPt ? fMaxNPt : fMaxCPt; }
157  Double_t PtEmc() const { return fPtEmc ; }
158  Double_t PtSub() const { return fPtSub ; }
159  Double_t PtSubVect() const { return fPtSubVect ; }
160  Int_t TrackAt(Int_t idx) const { return fTrackIDs.At(idx) ; }
161 
162  // Background subtraction
163  Double_t PtSub(Double_t rho, Bool_t save = kFALSE) ;
164  Double_t PtSubVect(Double_t rho, Bool_t save = kFALSE) ;
165  TLorentzVector SubtractRhoVect(Double_t rho, Bool_t save = kFALSE);
166 
167  // Jet constituents
168  AliVCluster *Cluster(Int_t idx) const;
169  AliVCluster *ClusterAt(Int_t idx, TClonesArray *ca) const;
170  Int_t ContainsCluster(AliVCluster* cluster, TClonesArray* clusters) const;
171  Int_t ContainsCluster(Int_t ic) const;
172  AliVCluster *GetLeadingCluster(TClonesArray *clusters = 0) const;
173  AliVParticle *Track(Int_t idx) const;
174  AliVParticle *TrackAt(Int_t idx, TClonesArray *ta) const;
175  Int_t ContainsTrack(AliVParticle* track, TClonesArray* tracks) const;
176  Int_t ContainsTrack(Int_t it) const;
177  AliVParticle *GetLeadingTrack(TClonesArray *tracks = 0) const;
178  Bool_t IsGhost() const { return fNn + fNch == 0; }
179 
184  const std::vector<PWG::JETFW::AliEmcalParticleJetConstituent> &GetParticleConstituents() const { return fParticleConstituents; }
185 
190  const std::vector<PWG::JETFW::AliEmcalClusterJetConstituent> &GetClusterConstituents() const { return fClusterConstituents; }
191 
197 
203 
210 
216  const PWG::JETFW::AliEmcalParticleJetConstituent *ParticleConstituentAt(unsigned int ipart) const;
217 
223 
229 
235  bool HasClusterConstituent(const AliVCluster *const clust) const;
236 
242  bool HasParticleConstituent(const AliVParticle *const part) const;
243 
244  // Fragmentation function
245  Double_t GetZ(const Double_t trkPx, const Double_t trkPy, const Double_t trkPz) const;
246  Double_t GetZ(const AliVParticle* trk ) const;
247  Double_t GetXi(const AliVParticle* trk ) const;
248  Double_t GetXi(const Double_t trkPx, const Double_t trkPy, const Double_t trkPz) const;
249 
250  // Other service methods
251  void GetMomentum(TLorentzVector &vec) const;
252  Double_t DeltaR(const AliVParticle* part) const;
253 
254  // Setters
255  void SetLabel(Int_t l) { fLabel = l; }
256  void SetArea(Double_t a) { fArea = a; }
257  void SetAreaEta(Double_t a) { fAreaEta = a; }
258  void SetAreaPhi(Double_t a) { fAreaPhi = TVector2::Phi_0_2pi(a); }
259  void SetAreaE(Double_t a) { fAreaE = a; }
260  void SetAreaEmc(Double_t a) { fAreaEmc = a; }
262  void SetMaxNeutralPt(Double32_t t) { fMaxNPt = t; }
263  void SetMaxChargedPt(Double32_t t) { fMaxCPt = t; }
264  void SetNEF(Double_t nef) { fNEF = nef; }
266  void SetNumberOfTracks(Int_t n) { fTrackIDs.Set(n); }
267  void SetNumberOfCharged(Int_t n) { fNch = n; }
268  void SetNumberOfNeutrals(Int_t n) { fNn = n; }
269  void SetMCPt(Double_t p) { fMCPt = p; }
270  void SetNEmc(Int_t n) { fNEmc = n; }
271  void SetPtEmc(Double_t pt) { fPtEmc = pt; }
272  void SetPtSub(Double_t ps) { fPtSub = ps; }
273  void SetPtSubVect(Double_t ps) { fPtSubVect = ps; }
274  void AddClusterAt(Int_t clus, Int_t idx){ fClusterIDs.AddAt(clus, idx); }
275  void AddTrackAt(Int_t track, Int_t idx) { fTrackIDs.AddAt(track, idx); }
276  void Clear(Option_t */*option*/="");
277 
285  void AddParticleConstituent(const AliVParticle *const part, Bool_t isFromEmbeddedEvent, UInt_t globalIndex);
286 
293 
303  void AddClusterConstituent(const AliVCluster *const clust, AliVCluster::VCluUserDefEnergy_t endef, Double_t *pvec, Bool_t isFromEmbeddedEvent, UInt_t globalIndex);
304 
311 
312  // Sorting methods
313  void SortConstituents();
314  std::vector<int> GetPtSortedTrackConstituentIndexes(TClonesArray *tracks) const;
315 
316  // Trigger
317  Bool_t IsTriggerJet(UInt_t trigger=AliVEvent::kEMCEJE) const { return (Bool_t)((fTriggers & trigger) != 0); }
318  void SetTrigger(UInt_t trigger) { fTriggers = trigger; }
319  void AddTrigger(UInt_t trigger) { fTriggers |= trigger; }
320 
321  // Matching
322  void ResetMatching();
327  AliEmcalJet* ClosestJet() const { return fClosestJets[0] ; }
329  AliEmcalJet* SecondClosestJet() const { return fClosestJets[1] ; }
331  AliEmcalJet* MatchedJet() const { return fMatched < 2 ? fClosestJets[fMatched] : 0; }
333 
334  // Jet tagging
336  void SetTagStatus(Int_t i) { fTagStatus = i ; }
337  AliEmcalJet* GetTaggedJet() const { return fTaggedJet ; }
338  Int_t GetTagStatus() const { return fTagStatus ; }
339 
340  // Ghosts
341  void AddGhost(const Double_t dPx, const Double_t dPy, const Double_t dPz, const Double_t dE);
342  Bool_t HasGhost() const { return fHasGhost; }
343  const std::vector<TLorentzVector> GetGhosts() const { return fGhosts ; }
344 
345  // Debug printouts
346  void Print(Option_t* /*opt*/ = "") const;
347  void PrintConstituents(TClonesArray* tracks = 0, TClonesArray* clusters = 0) const;
348 
349  //heavy-flavor jets
350  Int_t GetFlavour() const { return fFlavourTagging; }
351  void AddFlavourTag(Int_t tag) { fFlavourTagging |= tag; }
352  void SetFlavour(Int_t flavour) { fFlavourTagging = flavour; }
353  Bool_t TestFlavourTag(Int_t tag) const { return (Bool_t)((tag & fFlavourTagging) !=0); }
355  void AddFlavourTrack(AliVParticle* hftrack);
356  AliVParticle *GetFlavourTrack(Int_t i=0) const;
357  Double_t GetFlavourTrackZ(Int_t i=0) const;
358  AliVParticle *RemoveFlavourTrack(Int_t i=0) ;
359 
360  // Jet shape
364 
365  // Jet geometrical acceptance
368 
369  protected:
371  Double32_t fPt; //[0,0,12]
373  Double32_t fEta; //[-1,1,12]
375  Double32_t fPhi; //[0,6.3,12]
377  Double32_t fM; //[0,0,8]
379  Double32_t fNEF; //[0,1,8]
381  Double32_t fArea; //[0,0,12]
383  Double32_t fAreaEta; //[0,0,12]
385  Double32_t fAreaPhi; //[0,0,12]
387  Double32_t fAreaE; //[0,0,12]
389  Double32_t fAreaEmc; //[0,0,12]
393  Double32_t fMaxCPt; //[0,0,12]
396  Double32_t fMaxNPt; //[0,0,12]
397  Double32_t fMCPt;
400  Double32_t fPtEmc; //[0,0,12]
406  Double32_t fClosestJetsDist[2];
415 
417  std::vector<TLorentzVector> fGhosts;
418 
421 
422  std::vector<PWG::JETFW::AliEmcalParticleJetConstituent> fParticleConstituents;
423  std::vector<PWG::JETFW::AliEmcalClusterJetConstituent> fClusterConstituents;
424 
425  private:
430  struct sort_descend {
431  // first value of the pair is Pt and the second is entry index
432  bool operator () (const std::pair<Double_t, Int_t>& p1, const std::pair<Double_t, Int_t>& p2) { return p1.first > p2.first ; }
433  };
434 
436  ClassDef(AliEmcalJet,19);
438 };
439 
440 std::ostream &operator<<(std::ostream &in, const AliEmcalJet &jet);
441 #endif
TArrayI fTrackIDs
Array containing ids of track constituents.
Definition: AliEmcalJet.h:404
Double_t AreaEmc() const
Definition: AliEmcalJet.h:135
void SetMaxNeutralPt(Double32_t t)
Definition: AliEmcalJet.h:262
Int_t PdgCode() const
Definition: AliEmcalJet.h:125
void SetSecondClosestJet(AliEmcalJet *j, Double_t d)
Definition: AliEmcalJet.h:324
void SetTagStatus(Int_t i)
Definition: AliEmcalJet.h:336
Double_t fPtSubVect
! Background vector subtracted pt (not stored set from outside)
Definition: AliEmcalJet.h:412
Int_t NEmc() const
Definition: AliEmcalJet.h:152
const PWG::JETFW::AliEmcalParticleJetConstituent * ParticleConstituentAt(unsigned int ipart) const
Access to the -particle constituent.
Short_t Charge() const
Definition: AliEmcalJet.h:123
void AddFlavourTag(Int_t tag)
Definition: AliEmcalJet.h:351
Double_t Area() const
Definition: AliEmcalJet.h:130
void AddTrackAt(Int_t track, Int_t idx)
Definition: AliEmcalJet.h:275
Generic signal 2.
Definition: AliEmcalJet.h:88
Bool_t XvYvZv(Double_t x[3]) const
Definition: AliEmcalJet.h:115
double Double_t
Definition: External.C:58
AliEmcalJet * GetTaggedJet() const
Definition: AliEmcalJet.h:337
Double_t Xv() const
Definition: AliEmcalJet.h:112
Int_t fTagStatus
! Status of tagging -1: NA 0: not tagged 1: tagged
Definition: AliEmcalJet.h:410
Bool_t IsMC() const
Definition: AliEmcalJet.h:144
AliVParticle * GetLeadingTrack(TClonesArray *tracks=0) const
void SetTaggedJet(AliEmcalJet *j)
Definition: AliEmcalJet.h:335
Double_t GetXi(const AliVParticle *trk) const
Double_t MCPt() const
Definition: AliEmcalJet.h:153
Simple C structure to allow sorting in descending order.
Definition: AliEmcalJet.h:430
Double32_t fAreaE
Jet temporal area component.
Definition: AliEmcalJet.h:387
AliEmcalJet * ClosestJet() const
Definition: AliEmcalJet.h:327
Jet is tagged to contain a D* meson.
Definition: AliEmcalJet.h:85
Int_t GetTagStatus() const
Definition: AliEmcalJet.h:338
PHOS acceptance.
Definition: AliEmcalJet.h:75
Double_t ClosestJetDistance() const
Definition: AliEmcalJet.h:328
TArrayI fClusterIDs
Array containing ids of cluster constituents.
Definition: AliEmcalJet.h:403
DCal acceptance – spans entire rectangular region in eta-phi (including most of PHOS) ...
Definition: AliEmcalJet.h:71
Bool_t AxisInEmcal() const
Definition: AliEmcalJet.h:136
void ClearFlavourTracks()
Definition: AliEmcalJet.h:354
Double_t Eta() const
Definition: AliEmcalJet.h:121
Int_t fLabel
! Label to inclusive jet for constituent subtracted jet
Definition: AliEmcalJet.h:414
Double_t fPtSub
! Background subtracted pt (not stored set from outside)
Definition: AliEmcalJet.h:411
Int_t GetFlavour() const
Definition: AliEmcalJet.h:350
Double_t Py() const
Definition: AliEmcalJet.h:107
void Clear(Option_t *="")
Double_t Phi() const
Definition: AliEmcalJet.h:117
JetAcceptanceType
Bit definition for jet geometry acceptance. Cut implemented in AliJetContainer by comparing jet&#39;s bit...
Definition: AliEmcalJet.h:66
Bool_t fHasGhost
! Whether ghost particle are included within the constituents
Definition: AliEmcalJet.h:416
void AddClusterAt(Int_t clus, Int_t idx)
Definition: AliEmcalJet.h:274
std::vector< PWG::JETFW::AliEmcalParticleJetConstituent > fParticleConstituents
List of particle constituents.
Definition: AliEmcalJet.h:422
AliEmcalJet * MatchedJet() const
Definition: AliEmcalJet.h:331
Double_t AreaE() const
Definition: AliEmcalJet.h:134
Bool_t IsSortable() const
Definition: AliEmcalJet.h:145
Int_t GetLabel() const
Definition: AliEmcalJet.h:124
Int_t ClusterAt(Int_t idx) const
Definition: AliEmcalJet.h:137
void SetArea(Double_t a)
Definition: AliEmcalJet.h:256
AliEmcalJetShapeProperties * fJetShapeProperties
! Pointer to the jet shape properties
Definition: AliEmcalJet.h:419
AliVParticle * Track(Int_t idx) const
Full acceptance, i.e. no acceptance cut applied – left to user.
Definition: AliEmcalJet.h:77
Double_t FracEmcalArea() const
Definition: AliEmcalJet.h:141
const std::vector< PWG::JETFW::AliEmcalParticleJetConstituent > & GetParticleConstituents() const
Get container with particle (track / MC particle) constituents.
Definition: AliEmcalJet.h:184
Double_t MaxChargedPt() const
Definition: AliEmcalJet.h:147
AliEmcalJet & operator=(const AliEmcalJet &jet)
const PWG::JETFW::AliEmcalClusterJetConstituent * ClusterConstituentAt(unsigned int icl) const
Access to the -cluster constituent.
Double_t GetFlavourTrackZ(Int_t i=0) const
void SetAreaE(Double_t a)
Definition: AliEmcalJet.h:259
std::vector< PWG::JETFW::AliEmcalClusterJetConstituent > fClusterConstituents
List of cluster constituents.
Definition: AliEmcalJet.h:423
void SetMCPt(Double_t p)
Definition: AliEmcalJet.h:269
Double_t E() const
Definition: AliEmcalJet.h:119
const PWG::JETFW::AliEmcalParticleJetConstituent * GetLeadingParticleConstituent() const
Get the leading particle constituent.
void SetPtSubVect(Double_t ps)
Definition: AliEmcalJet.h:273
UShort_t Nch() const
Definition: AliEmcalJet.h:150
AliVParticle * GetFlavourTrack(Int_t i=0) const
void SetLabel(Int_t l)
Definition: AliEmcalJet.h:255
Double_t AreaPt() const
Definition: AliEmcalJet.h:131
UShort_t GetNumberOfConstituents() const
Definition: AliEmcalJet.h:140
AliVParticle * RemoveFlavourTrack(Int_t i=0)
Int_t ContainsCluster(AliVCluster *cluster, TClonesArray *clusters) const
void SetMatchedToSecondClosest(UShort_t m)
Definition: AliEmcalJet.h:326
Int_t TrackAt(Int_t idx) const
Definition: AliEmcalJet.h:160
UShort_t GetNumberOfTracks() const
Definition: AliEmcalJet.h:139
AliEmcalJetShapeProperties * GetShapeProperties()
Definition: AliEmcalJet.h:362
Double_t Px() const
Definition: AliEmcalJet.h:106
void AddParticleConstituent(const AliVParticle *const part, Bool_t isFromEmbeddedEvent, UInt_t globalIndex)
Add new particle (track / mc particle) constituent to the given jet Note: this will append the consti...
std::ostream & Print(std::ostream &in) const
void ResetMatching()
void GetMomentum(TLorentzVector &vec) const
Double32_t fMCPt
Pt from MC particles contributing to the jet.
Definition: AliEmcalJet.h:397
Double32_t fMaxNPt
Pt of maximum neutral constituent.
Definition: AliEmcalJet.h:396
TPC fiducial acceptance (each eta edge narrowed by jet R)
Definition: AliEmcalJet.h:68
void SetJetAcceptanceType(UInt_t type)
Definition: AliEmcalJet.h:366
Double32_t fM
Jet mass.
Definition: AliEmcalJet.h:377
Double_t AreaEta() const
Definition: AliEmcalJet.h:132
Bool_t IsInEmcal() const
Definition: AliEmcalJet.h:143
int Int_t
Definition: External.C:63
UShort_t GetNumberOfClusters() const
Definition: AliEmcalJet.h:138
friend std::ostream & operator<<(std::ostream &in, const AliEmcalJet &myjet)
TPC acceptance.
Definition: AliEmcalJet.h:67
unsigned int UInt_t
Definition: External.C:33
UShort_t GetMatchingType() const
Definition: AliEmcalJet.h:332
Double32_t fAreaPhi
Jet phi area.
Definition: AliEmcalJet.h:385
Double_t Yv() const
Definition: AliEmcalJet.h:113
UInt_t GetJetAcceptanceType() const
Definition: AliEmcalJet.h:367
PHOS fiducial acceptance (each eta, phi edge narrowed by jet R)
Definition: AliEmcalJet.h:76
void AddFlavourTrack(AliVParticle *hftrack)
EFlavourTag
Bit definition for the flavor tagging.
Definition: AliEmcalJet.h:84
void AddClusterConstituent(const AliVCluster *const clust, AliVCluster::VCluUserDefEnergy_t endef, Double_t *pvec, Bool_t isFromEmbeddedEvent, UInt_t globalIndex)
Add new cluster constituent to the given jet Note: this will append the constituent. No sorting according to particle is done.
EMCal acceptance.
Definition: AliEmcalJet.h:69
UShort_t fMatched
! 0 or 1 if it is matched with one of the closest jets; 2 if it is not matched
Definition: AliEmcalJet.h:407
std::vector< TLorentzVector > fGhosts
! Vector containing the ghost particles
Definition: AliEmcalJet.h:417
Double_t OneOverPt() const
Definition: AliEmcalJet.h:116
TString toString() const
Double_t PtSubVect() const
Definition: AliEmcalJet.h:159
Bool_t IsGhost() const
Definition: AliEmcalJet.h:178
Double32_t fArea
Jet transverse area.
Definition: AliEmcalJet.h:381
Double_t Theta() const
Definition: AliEmcalJet.h:118
Double_t MaxTrackPt() const
Definition: AliEmcalJet.h:155
AliEmcalJet * fTaggedJet
! Jet tagged to this jet
Definition: AliEmcalJet.h:409
Int_t ContainsTrack(AliVParticle *track, TClonesArray *tracks) const
Double_t Phi_0_2pi() const
Definition: AliEmcalJet.h:129
TLorentzVector SubtractRhoVect(Double_t rho, Bool_t save=kFALSE)
Generic background 1.
Definition: AliEmcalJet.h:89
virtual ~AliEmcalJet()
Double_t PtSub() const
Definition: AliEmcalJet.h:158
int GetNumberOfClusterConstituents() const
Get the number of cluster constituents.
Definition: AliEmcalJet.h:202
bool HasParticleConstituent(const AliVParticle *const part) const
Checks whether a given particle is a jet constituent.
Double_t Y() const
Definition: AliEmcalJet.h:122
AliVCluster * GetLeadingCluster(TClonesArray *clusters=0) const
Double32_t fAreaEta
Jet eta area.
Definition: AliEmcalJet.h:383
void SetPtSub(Double_t ps)
Definition: AliEmcalJet.h:272
void SetMaxChargedPt(Double32_t t)
Definition: AliEmcalJet.h:263
void SetNEF(Double_t nef)
Definition: AliEmcalJet.h:264
const Double_t * PID() const
Definition: AliEmcalJet.h:126
Double_t MaxNeutralPt() const
Definition: AliEmcalJet.h:146
void AddGhost(const Double_t dPx, const Double_t dPy, const Double_t dPz, const Double_t dE)
Double_t GetZ(const Double_t trkPx, const Double_t trkPy, const Double_t trkPz) const
void SetNEmc(Int_t n)
Definition: AliEmcalJet.h:270
void PrintConstituents(TClonesArray *tracks=0, TClonesArray *clusters=0) const
Double_t PtEmc() const
Definition: AliEmcalJet.h:157
Int_t fFlavourTagging
Tag jet with a flavor (use bits defined in enum EFlavourTag)
Definition: AliEmcalJet.h:391
AliEmcalJet * SecondClosestJet() const
Definition: AliEmcalJet.h:329
Implementation of a jet constituent for constituent clusters.
Double_t DeltaR(const AliVParticle *part) const
UShort_t N() const
Definition: AliEmcalJet.h:151
short Short_t
Definition: External.C:23
bool operator()(const std::pair< Double_t, Int_t > &p1, const std::pair< Double_t, Int_t > &p2)
Definition: AliEmcalJet.h:432
Double_t AreaPhi() const
Definition: AliEmcalJet.h:133
UInt_t fJetAcceptanceType
! Jet acceptance type (stored bitwise)
Definition: AliEmcalJet.h:420
Double_t Pt() const
Definition: AliEmcalJet.h:109
Double32_t fPtEmc
Pt in EMCAL acceptance.
Definition: AliEmcalJet.h:401
Bool_t HasGhost() const
Definition: AliEmcalJet.h:342
void SetNumberOfCharged(Int_t n)
Definition: AliEmcalJet.h:267
This class contains the derivative subtraction operators for jet shapes.
void SetAreaEta(Double_t a)
Definition: AliEmcalJet.h:257
const PWG::JETFW::AliEmcalClusterJetConstituent * GetLeadingClusterConstituent() const
Get the leading cluster constituent.
Jet is tagged to contain a D0 meson.
Definition: AliEmcalJet.h:86
Double32_t fPt
Jet transverse momentum.
Definition: AliEmcalJet.h:371
const std::vector< PWG::JETFW::AliEmcalClusterJetConstituent > & GetClusterConstituents() const
Get container with cluster constituents.
Definition: AliEmcalJet.h:190
Double_t P() const
Definition: AliEmcalJet.h:110
Bool_t IsInsideEmcal() const
Definition: AliEmcalJet.h:142
void SetClosestJet(AliEmcalJet *j, Double_t d)
Definition: AliEmcalJet.h:323
Double32_t fEta
Jet pseudo-rapidity.
Definition: AliEmcalJet.h:373
Double32_t fMaxCPt
Pt of maximum charged constituent.
Definition: AliEmcalJet.h:394
Double_t SecondClosestJetDistance() const
Definition: AliEmcalJet.h:330
Bool_t TestFlavourTag(Int_t tag) const
Definition: AliEmcalJet.h:353
Double32_t fAreaEmc
Area on EMCAL surface (determined by ghosts in EMCal acceptance)
Definition: AliEmcalJet.h:389
UShort_t Nn() const
Definition: AliEmcalJet.h:149
AliVCluster * Cluster(Int_t idx) const
int GetNumberOfParticleConstituents() const
Get the number of particle constituents assigned to the given jet.
Definition: AliEmcalJet.h:196
Double_t MaxClusterPt() const
Definition: AliEmcalJet.h:154
Represent a jet reconstructed using the EMCal jet framework.
Definition: AliEmcalJet.h:51
Int_t fNn
Number of neutral constituents.
Definition: AliEmcalJet.h:398
void SetNumberOfTracks(Int_t n)
Definition: AliEmcalJet.h:266
Bool_t PxPyPz(Double_t p[3]) const
Definition: AliEmcalJet.h:111
unsigned short UShort_t
Definition: External.C:28
Double_t Pz() const
Definition: AliEmcalJet.h:108
const std::vector< TLorentzVector > GetGhosts() const
Definition: AliEmcalJet.h:343
void SetTrigger(UInt_t trigger)
Definition: AliEmcalJet.h:318
void SetFlavour(Int_t flavour)
Definition: AliEmcalJet.h:352
const char Option_t
Definition: External.C:48
UInt_t fTriggers
! Triggers that the jet might have fired (AliVEvent::EOfflineTriggerTypes)
Definition: AliEmcalJet.h:413
Bool_t fAxisInEmcal
Whether the jet axis is inside the EMCAL acceptance.
Definition: AliEmcalJet.h:390
void SortConstituents()
Int_t fNEmc
Number of constituents in EMCAL acceptance.
Definition: AliEmcalJet.h:402
void CreateShapeProperties()
Definition: AliEmcalJet.h:363
Generic signal 1.
Definition: AliEmcalJet.h:87
Generic background 2.
Definition: AliEmcalJet.h:90
Double32_t fPhi
Jet axis azimuthal angle.
Definition: AliEmcalJet.h:375
DCal fiducial acceptance (each eta, phi edge narrowed by jet R)
Definition: AliEmcalJet.h:72
AliEmcalJet * fClosestJets[2]
! If this is MC it contains the two closest detector level jets in order of distance and viceversa ...
Definition: AliEmcalJet.h:405
bool Bool_t
Definition: External.C:53
AliEmcalJetShapeProperties * GetShapeProperties() const
Definition: AliEmcalJet.h:361
void AddTrigger(UInt_t trigger)
Definition: AliEmcalJet.h:319
DCal acceptance – spans ONLY DCal (no PHOS or gap)
Definition: AliEmcalJet.h:73
Double_t NEF() const
Definition: AliEmcalJet.h:148
void SetAreaPhi(Double_t a)
Definition: AliEmcalJet.h:258
void SetPtEmc(Double_t pt)
Definition: AliEmcalJet.h:271
Bool_t IsTriggerJet(UInt_t trigger=AliVEvent::kEMCEJE) const
Definition: AliEmcalJet.h:317
bool HasClusterConstituent(const AliVCluster *const clust) const
Checks whether a given cluster is a constituent of the jet.
Double32_t fNEF
Jet Neutral Energy Fraction.
Definition: AliEmcalJet.h:379
Double32_t fClosestJetsDist[2]
! Distance from the two closest jets
Definition: AliEmcalJet.h:406
DCal fiducial acceptance (each eta, phi edge narrowed by jet R)
Definition: AliEmcalJet.h:74
UShort_t fMatchingType
! Matching type
Definition: AliEmcalJet.h:408
Double_t M() const
Definition: AliEmcalJet.h:120
EMCal fiducial acceptance (each eta, phi edge narrowed by jet R)
Definition: AliEmcalJet.h:70
Int_t Compare(const TObject *obj) const
Double_t MaxPartPt() const
Definition: AliEmcalJet.h:156
Generic background 3.
Definition: AliEmcalJet.h:91
void SetMatchedToClosest(UShort_t m)
Definition: AliEmcalJet.h:325
void SetNumberOfClusters(Int_t n)
Definition: AliEmcalJet.h:265
std::vector< int > GetPtSortedTrackConstituentIndexes(TClonesArray *tracks) const
Double_t Zv() const
Definition: AliEmcalJet.h:114
void SetAxisInEmcal(Bool_t b)
Definition: AliEmcalJet.h:261
TObjArray * fFlavourTracks
Definition: AliEmcalJet.h:392
void SetAreaEmc(Double_t a)
Definition: AliEmcalJet.h:260
void SetNumberOfNeutrals(Int_t n)
Definition: AliEmcalJet.h:268