AliPhysics  19b3b9d (19b3b9d)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 
23 
44 class AliEmcalJet : public AliVParticle
45 {
46  public:
47 
60  kTPC = 1<<0,
61  kTPCfid = 1<<1,
62  kEMCAL = 1<<2,
63  kEMCALfid = 1<<3,
64  kDCAL = 1<<4,
65  kDCALfid = 1<<5,
66  kDCALonly = 1<<6,
67  kDCALonlyfid = 1<<7,
68  kPHOS = 1<<8,
69  kPHOSfid = 1<<9,
70  kUser = 1<<10
71  };
72 
77  enum EFlavourTag {
78  kDStar = 1<<0,
79  kD0 = 1<<1,
80  kSig1 = 1<<2,
81  kSig2 = 1<<3,
82  kBckgrd1 = 1<<4,
83  kBckgrd2 = 1<<5,
84  kBckgrd3 = 1<<6
85  };
86 
87  AliEmcalJet();
90  AliEmcalJet(const AliEmcalJet &jet);
91  AliEmcalJet& operator=(const AliEmcalJet &jet);
92  virtual ~AliEmcalJet();
93  friend std::ostream &operator<<(std::ostream &in, const AliEmcalJet &jet);
94  Int_t Compare(const TObject* obj) const;
95  std::ostream &Print(std::ostream &in) const;
96  TString toString() const;
97 
98  // Implementation of AliVParticle interface
99  Double_t Px() const { return fPt*TMath::Cos(fPhi) ; }
100  Double_t Py() const { return fPt*TMath::Sin(fPhi) ; }
101  Double_t Pz() const { return fPt*TMath::SinH(fEta); }
102  Double_t Pt() const { return fPt ; }
103  Double_t P() const { return fPt*TMath::CosH(fEta); }
104  Bool_t PxPyPz(Double_t p[3]) const { p[0]=Px();p[1]=Py();p[2]=Pz(); return kTRUE; }
105  Double_t Xv() const { return 0.; }
106  Double_t Yv() const { return 0.; }
107  Double_t Zv() const { return 0.; }
108  Bool_t XvYvZv(Double_t x[3]) const { x[0]=0;x[1]=0;x[2]=0 ; return kTRUE; }
109  Double_t OneOverPt() const { return 1./fPt ; }
110  Double_t Phi() const { return fPhi ; }
111  Double_t Theta() const { return 2*TMath::ATan(TMath::Exp(-fEta)) ; }
112  Double_t E() const { Double_t p=P(); return TMath::Sqrt(fM*fM+p*p); }
113  Double_t M() const { return fM ; }
114  Double_t Eta() const { return fEta ; }
115  Double_t Y() const { Double_t e = E(); Double_t pz = Pz(); return 0.5*TMath::Log((e+pz)/(e-pz)); }
116  Short_t Charge() const { return 0 ; }
117  Int_t GetLabel() const { return fLabel ; }
118  Int_t PdgCode() const { return 0; }
119  const Double_t *PID() const { return 0; }
120 
121  // Other kinematic and jet properties
122  Double_t Phi_0_2pi() const { return TVector2::Phi_0_2pi(fPhi); }
123  Double_t Area() const { return fArea ; }
124  Double_t AreaPt() const { return fArea ; }
125  Double_t AreaEta() const { return fAreaEta ; }
126  Double_t AreaPhi() const { return fAreaPhi ; }
127  Double_t AreaE() const { return fAreaE ; }
128  Double_t AreaEmc() const { return fAreaEmc ; }
129  Bool_t AxisInEmcal() const { return fAxisInEmcal ; }
130  Int_t ClusterAt(Int_t idx) const { return fClusterIDs.At(idx) ; }
131  UShort_t GetNumberOfClusters() const { return fClusterIDs.GetSize() ; }
132  UShort_t GetNumberOfTracks() const { return fTrackIDs.GetSize() ; }
134  Double_t FracEmcalArea() const { return fAreaEmc/fArea ; }
135  Bool_t IsInsideEmcal() const { return (fAreaEmc/fArea>0.999) ; }
136  Bool_t IsInEmcal() const { return (Bool_t)(fAreaEmc > 0) ; }
137  Bool_t IsMC() const { return (Bool_t)(MCPt() > 0) ; }
138  Bool_t IsSortable() const { return kTRUE ; }
139  Double_t MaxNeutralPt() const { return fMaxNPt ; }
140  Double_t MaxChargedPt() const { return fMaxCPt ; }
141  Double_t NEF() const { return fNEF ; }
142  UShort_t Nn() const { return fNn ; }
143  UShort_t Nch() const { return fNch ; }
144  UShort_t N() const { return Nch()+Nn() ; }
145  Int_t NEmc() const { return fNEmc ; }
146  Double_t MCPt() const { return fMCPt ; }
147  Double_t MaxClusterPt() const { return MaxNeutralPt() ; }
148  Double_t MaxTrackPt() const { return MaxChargedPt() ; }
149  Double_t MaxPartPt() const { return fMaxCPt < fMaxNPt ? fMaxNPt : fMaxCPt; }
150  Double_t PtEmc() const { return fPtEmc ; }
151  Double_t PtSub() const { return fPtSub ; }
152  Double_t PtSubVect() const { return fPtSubVect ; }
153  Int_t TrackAt(Int_t idx) const { return fTrackIDs.At(idx) ; }
154 
155  // Background subtraction
156  Double_t PtSub(Double_t rho, Bool_t save = kFALSE) ;
157  Double_t PtSubVect(Double_t rho, Bool_t save = kFALSE) ;
158  TLorentzVector SubtractRhoVect(Double_t rho, Bool_t save = kFALSE);
159 
160  // Jet constituents
161  AliVCluster *Cluster(Int_t idx) const;
162  AliVCluster *ClusterAt(Int_t idx, TClonesArray *ca) const;
163  Int_t ContainsCluster(AliVCluster* cluster, TClonesArray* clusters) const;
164  Int_t ContainsCluster(Int_t ic) const;
165  AliVCluster *GetLeadingCluster(TClonesArray *clusters) const;
166  AliVParticle *Track(Int_t idx) const;
167  AliVParticle *TrackAt(Int_t idx, TClonesArray *ta) const;
168  Int_t ContainsTrack(AliVParticle* track, TClonesArray* tracks) const;
169  Int_t ContainsTrack(Int_t it) const;
170  AliVParticle *GetLeadingTrack(TClonesArray *tracks) const;
171  Bool_t IsGhost() const { return fNn + fNch == 0; }
172 
173  // Fragmentation function
174  Double_t GetZ(const Double_t trkPx, const Double_t trkPy, const Double_t trkPz) const;
175  Double_t GetZ(const AliVParticle* trk ) const;
176  Double_t GetXi(const AliVParticle* trk ) const;
177  Double_t GetXi(const Double_t trkPx, const Double_t trkPy, const Double_t trkPz) const;
178 
179  // Other service methods
180  void GetMomentum(TLorentzVector &vec) const;
181  Double_t DeltaR(const AliVParticle* part) const;
182 
183  // Setters
184  void SetLabel(Int_t l) { fLabel = l; }
185  void SetArea(Double_t a) { fArea = a; }
186  void SetAreaEta(Double_t a) { fAreaEta = a; }
187  void SetAreaPhi(Double_t a) { fAreaPhi = TVector2::Phi_0_2pi(a); }
188  void SetAreaE(Double_t a) { fAreaE = a; }
189  void SetAreaEmc(Double_t a) { fAreaEmc = a; }
191  void SetMaxNeutralPt(Double32_t t) { fMaxNPt = t; }
192  void SetMaxChargedPt(Double32_t t) { fMaxCPt = t; }
193  void SetNEF(Double_t nef) { fNEF = nef; }
195  void SetNumberOfTracks(Int_t n) { fTrackIDs.Set(n); }
196  void SetNumberOfCharged(Int_t n) { fNch = n; }
197  void SetNumberOfNeutrals(Int_t n) { fNn = n; }
198  void SetMCPt(Double_t p) { fMCPt = p; }
199  void SetNEmc(Int_t n) { fNEmc = n; }
200  void SetPtEmc(Double_t pt) { fPtEmc = pt; }
201  void SetPtSub(Double_t ps) { fPtSub = ps; }
202  void SetPtSubVect(Double_t ps) { fPtSubVect = ps; }
203  void AddClusterAt(Int_t clus, Int_t idx){ fClusterIDs.AddAt(clus, idx); }
204  void AddTrackAt(Int_t track, Int_t idx) { fTrackIDs.AddAt(track, idx); }
205  void Clear(Option_t */*option*/="");
206 
207  // Sorting methods
208  void SortConstituents();
209  std::vector<int> GetPtSortedTrackConstituentIndexes(TClonesArray *tracks) const;
210 
211  // Trigger
212  Bool_t IsTriggerJet(UInt_t trigger=AliVEvent::kEMCEJE) const { return (Bool_t)((fTriggers & trigger) != 0); }
213  void SetTrigger(UInt_t trigger) { fTriggers = trigger; }
214  void AddTrigger(UInt_t trigger) { fTriggers |= trigger; }
215 
216  // Matching
217  void ResetMatching();
222  AliEmcalJet* ClosestJet() const { return fClosestJets[0] ; }
224  AliEmcalJet* SecondClosestJet() const { return fClosestJets[1] ; }
226  AliEmcalJet* MatchedJet() const { return fMatched < 2 ? fClosestJets[fMatched] : 0; }
228 
229  // Jet tagging
231  void SetTagStatus(Int_t i) { fTagStatus = i ; }
232  AliEmcalJet* GetTaggedJet() const { return fTaggedJet ; }
233  Int_t GetTagStatus() const { return fTagStatus ; }
234 
235  // Ghosts
236  void AddGhost(const Double_t dPx, const Double_t dPy, const Double_t dPz, const Double_t dE);
237  Bool_t HasGhost() const { return fHasGhost; }
238  const std::vector<TLorentzVector> GetGhosts() const { return fGhosts ; }
239 
240  // Debug printouts
241  void Print(Option_t* /*opt*/ = "") const;
242  void PrintConstituents(TClonesArray* tracks, TClonesArray* clusters) const;
243 
244  //heavy-flavor jets
245  Int_t GetFlavour() const { return fFlavourTagging; }
246  void AddFlavourTag(Int_t tag) { fFlavourTagging |= tag; }
247  void SetFlavour(Int_t flavour) { fFlavourTagging = flavour; }
248  Bool_t TestFlavourTag(Int_t tag) const { return (Bool_t)((tag & fFlavourTagging) !=0); }
250  void AddFlavourTrack(AliVParticle* hftrack);
251  AliVParticle *GetFlavourTrack(Int_t i=0) const;
252  Double_t GetFlavourTrackZ(Int_t i=0) const;
253  AliVParticle *RemoveFlavourTrack(Int_t i=0) ;
254 
255  // Jet shape
259 
260  // Jet geometrical acceptance
263 
264  protected:
266  Double32_t fPt; //[0,0,12]
268  Double32_t fEta; //[-1,1,12]
270  Double32_t fPhi; //[0,6.3,12]
272  Double32_t fM; //[0,0,8]
274  Double32_t fNEF; //[0,1,8]
276  Double32_t fArea; //[0,0,12]
278  Double32_t fAreaEta; //[0,0,12]
280  Double32_t fAreaPhi; //[0,0,12]
282  Double32_t fAreaE; //[0,0,12]
284  Double32_t fAreaEmc; //[0,0,12]
288  Double32_t fMaxCPt; //[0,0,12]
291  Double32_t fMaxNPt; //[0,0,12]
292  Double32_t fMCPt;
295  Double32_t fPtEmc; //[0,0,12]
301  Double32_t fClosestJetsDist[2];
310 
312  std::vector<TLorentzVector> fGhosts;
313 
316 
317  private:
322  struct sort_descend {
323  // first value of the pair is Pt and the second is entry index
324  bool operator () (const std::pair<Double_t, Int_t>& p1, const std::pair<Double_t, Int_t>& p2) { return p1.first > p2.first ; }
325  };
326 
328  ClassDef(AliEmcalJet,19);
330 };
331 
332 std::ostream &operator<<(std::ostream &in, const AliEmcalJet &jet);
333 #endif
TArrayI fTrackIDs
Array containing ids of track constituents.
Definition: AliEmcalJet.h:299
Double_t AreaEmc() const
Definition: AliEmcalJet.h:128
void SetMaxNeutralPt(Double32_t t)
Definition: AliEmcalJet.h:191
Int_t PdgCode() const
Definition: AliEmcalJet.h:118
void SetSecondClosestJet(AliEmcalJet *j, Double_t d)
Definition: AliEmcalJet.h:219
void SetTagStatus(Int_t i)
Definition: AliEmcalJet.h:231
Double_t fPtSubVect
! Background vector subtracted pt (not stored set from outside)
Definition: AliEmcalJet.h:307
Int_t NEmc() const
Definition: AliEmcalJet.h:145
Short_t Charge() const
Definition: AliEmcalJet.h:116
void AddFlavourTag(Int_t tag)
Definition: AliEmcalJet.h:246
Double_t Area() const
Definition: AliEmcalJet.h:123
void AddTrackAt(Int_t track, Int_t idx)
Definition: AliEmcalJet.h:204
Generic signal 2.
Definition: AliEmcalJet.h:81
Bool_t XvYvZv(Double_t x[3]) const
Definition: AliEmcalJet.h:108
double Double_t
Definition: External.C:58
AliEmcalJet * GetTaggedJet() const
Definition: AliEmcalJet.h:232
Double_t Xv() const
Definition: AliEmcalJet.h:105
Int_t fTagStatus
! Status of tagging -1: NA 0: not tagged 1: tagged
Definition: AliEmcalJet.h:305
Bool_t IsMC() const
Definition: AliEmcalJet.h:137
void SetTaggedJet(AliEmcalJet *j)
Definition: AliEmcalJet.h:230
Double_t GetXi(const AliVParticle *trk) const
Double_t MCPt() const
Definition: AliEmcalJet.h:146
Simple C structure to allow sorting in descending order.
Definition: AliEmcalJet.h:322
Double32_t fAreaE
Jet temporal area component.
Definition: AliEmcalJet.h:282
AliEmcalJet * ClosestJet() const
Definition: AliEmcalJet.h:222
Jet is tagged to contain a D* meson.
Definition: AliEmcalJet.h:78
Int_t GetTagStatus() const
Definition: AliEmcalJet.h:233
PHOS acceptance.
Definition: AliEmcalJet.h:68
Double_t ClosestJetDistance() const
Definition: AliEmcalJet.h:223
TArrayI fClusterIDs
Array containing ids of cluster constituents.
Definition: AliEmcalJet.h:298
DCal acceptance – spans entire rectangular region in eta-phi (including most of PHOS) ...
Definition: AliEmcalJet.h:64
Bool_t AxisInEmcal() const
Definition: AliEmcalJet.h:129
void ClearFlavourTracks()
Definition: AliEmcalJet.h:249
Double_t Eta() const
Definition: AliEmcalJet.h:114
Int_t fLabel
! Label to inclusive jet for constituent subtracted jet
Definition: AliEmcalJet.h:309
Double_t fPtSub
! Background subtracted pt (not stored set from outside)
Definition: AliEmcalJet.h:306
Int_t GetFlavour() const
Definition: AliEmcalJet.h:245
Double_t Py() const
Definition: AliEmcalJet.h:100
void Clear(Option_t *="")
Double_t Phi() const
Definition: AliEmcalJet.h:110
JetAcceptanceType
Bit definition for jet geometry acceptance. Cut implemented in AliJetContainer by comparing jet's bit...
Definition: AliEmcalJet.h:59
Bool_t fHasGhost
! Whether ghost particle are included within the constituents
Definition: AliEmcalJet.h:311
void AddClusterAt(Int_t clus, Int_t idx)
Definition: AliEmcalJet.h:203
AliEmcalJet * MatchedJet() const
Definition: AliEmcalJet.h:226
Double_t AreaE() const
Definition: AliEmcalJet.h:127
Bool_t IsSortable() const
Definition: AliEmcalJet.h:138
Int_t GetLabel() const
Definition: AliEmcalJet.h:117
Int_t ClusterAt(Int_t idx) const
Definition: AliEmcalJet.h:130
void SetArea(Double_t a)
Definition: AliEmcalJet.h:185
AliEmcalJetShapeProperties * fJetShapeProperties
! Pointer to the jet shape properties
Definition: AliEmcalJet.h:314
AliVParticle * Track(Int_t idx) const
Full acceptance, i.e. no acceptance cut applied – left to user.
Definition: AliEmcalJet.h:70
Double_t FracEmcalArea() const
Definition: AliEmcalJet.h:134
Double_t MaxChargedPt() const
Definition: AliEmcalJet.h:140
AliEmcalJet & operator=(const AliEmcalJet &jet)
Double_t GetFlavourTrackZ(Int_t i=0) const
void SetAreaE(Double_t a)
Definition: AliEmcalJet.h:188
friend std::ostream & operator<<(std::ostream &in, const AliEmcalJet &jet)
void SetMCPt(Double_t p)
Definition: AliEmcalJet.h:198
Double_t E() const
Definition: AliEmcalJet.h:112
void SetPtSubVect(Double_t ps)
Definition: AliEmcalJet.h:202
UShort_t Nch() const
Definition: AliEmcalJet.h:143
AliVParticle * GetFlavourTrack(Int_t i=0) const
void SetLabel(Int_t l)
Definition: AliEmcalJet.h:184
Double_t AreaPt() const
Definition: AliEmcalJet.h:124
UShort_t GetNumberOfConstituents() const
Definition: AliEmcalJet.h:133
AliVParticle * RemoveFlavourTrack(Int_t i=0)
Int_t ContainsCluster(AliVCluster *cluster, TClonesArray *clusters) const
void SetMatchedToSecondClosest(UShort_t m)
Definition: AliEmcalJet.h:221
Int_t TrackAt(Int_t idx) const
Definition: AliEmcalJet.h:153
UShort_t GetNumberOfTracks() const
Definition: AliEmcalJet.h:132
AliEmcalJetShapeProperties * GetShapeProperties()
Definition: AliEmcalJet.h:257
Double_t Px() const
Definition: AliEmcalJet.h:99
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:292
Double32_t fMaxNPt
Pt of maximum neutral constituent.
Definition: AliEmcalJet.h:291
TPC fiducial acceptance (each eta edge narrowed by jet R)
Definition: AliEmcalJet.h:61
void SetJetAcceptanceType(UInt_t type)
Definition: AliEmcalJet.h:261
Double32_t fM
Jet mass.
Definition: AliEmcalJet.h:272
Double_t AreaEta() const
Definition: AliEmcalJet.h:125
Bool_t IsInEmcal() const
Definition: AliEmcalJet.h:136
int Int_t
Definition: External.C:63
UShort_t GetNumberOfClusters() const
Definition: AliEmcalJet.h:131
AliVParticle * GetLeadingTrack(TClonesArray *tracks) const
TPC acceptance.
Definition: AliEmcalJet.h:60
unsigned int UInt_t
Definition: External.C:33
UShort_t GetMatchingType() const
Definition: AliEmcalJet.h:227
Double32_t fAreaPhi
Jet phi area.
Definition: AliEmcalJet.h:280
Double_t Yv() const
Definition: AliEmcalJet.h:106
UInt_t GetJetAcceptanceType() const
Definition: AliEmcalJet.h:262
PHOS fiducial acceptance (each eta, phi edge narrowed by jet R)
Definition: AliEmcalJet.h:69
void AddFlavourTrack(AliVParticle *hftrack)
EFlavourTag
Bit definition for the flavor tagging.
Definition: AliEmcalJet.h:77
EMCal acceptance.
Definition: AliEmcalJet.h:62
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:302
std::vector< TLorentzVector > fGhosts
! Vector containing the ghost particles
Definition: AliEmcalJet.h:312
Double_t OneOverPt() const
Definition: AliEmcalJet.h:109
TString toString() const
Double_t PtSubVect() const
Definition: AliEmcalJet.h:152
Bool_t IsGhost() const
Definition: AliEmcalJet.h:171
Double32_t fArea
Jet transverse area.
Definition: AliEmcalJet.h:276
Double_t Theta() const
Definition: AliEmcalJet.h:111
Double_t MaxTrackPt() const
Definition: AliEmcalJet.h:148
AliVCluster * GetLeadingCluster(TClonesArray *clusters) const
AliEmcalJet * fTaggedJet
! Jet tagged to this jet
Definition: AliEmcalJet.h:304
Int_t ContainsTrack(AliVParticle *track, TClonesArray *tracks) const
Double_t Phi_0_2pi() const
Definition: AliEmcalJet.h:122
TLorentzVector SubtractRhoVect(Double_t rho, Bool_t save=kFALSE)
Generic background 1.
Definition: AliEmcalJet.h:82
virtual ~AliEmcalJet()
Double_t PtSub() const
Definition: AliEmcalJet.h:151
Double_t Y() const
Definition: AliEmcalJet.h:115
Double32_t fAreaEta
Jet eta area.
Definition: AliEmcalJet.h:278
void SetPtSub(Double_t ps)
Definition: AliEmcalJet.h:201
void SetMaxChargedPt(Double32_t t)
Definition: AliEmcalJet.h:192
void SetNEF(Double_t nef)
Definition: AliEmcalJet.h:193
const Double_t * PID() const
Definition: AliEmcalJet.h:119
Double_t MaxNeutralPt() const
Definition: AliEmcalJet.h:139
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:199
Double_t PtEmc() const
Definition: AliEmcalJet.h:150
Int_t fFlavourTagging
Tag jet with a flavor (use bits defined in enum EFlavourTag)
Definition: AliEmcalJet.h:286
AliEmcalJet * SecondClosestJet() const
Definition: AliEmcalJet.h:224
Double_t DeltaR(const AliVParticle *part) const
UShort_t N() const
Definition: AliEmcalJet.h:144
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:324
Double_t AreaPhi() const
Definition: AliEmcalJet.h:126
UInt_t fJetAcceptanceType
! Jet acceptance type (stored bitwise)
Definition: AliEmcalJet.h:315
Double_t Pt() const
Definition: AliEmcalJet.h:102
Double32_t fPtEmc
Pt in EMCAL acceptance.
Definition: AliEmcalJet.h:296
Bool_t HasGhost() const
Definition: AliEmcalJet.h:237
void SetNumberOfCharged(Int_t n)
Definition: AliEmcalJet.h:196
This class contains the derivative subtraction operators for jet shapes.
void SetAreaEta(Double_t a)
Definition: AliEmcalJet.h:186
Jet is tagged to contain a D0 meson.
Definition: AliEmcalJet.h:79
Double32_t fPt
Jet transverse momentum.
Definition: AliEmcalJet.h:266
Double_t P() const
Definition: AliEmcalJet.h:103
std::ostream & operator<<(std::ostream &in, const AliEmcalJet &jet)
Bool_t IsInsideEmcal() const
Definition: AliEmcalJet.h:135
void SetClosestJet(AliEmcalJet *j, Double_t d)
Definition: AliEmcalJet.h:218
Double32_t fEta
Jet pseudo-rapidity.
Definition: AliEmcalJet.h:268
Double32_t fMaxCPt
Pt of maximum charged constituent.
Definition: AliEmcalJet.h:289
Double_t SecondClosestJetDistance() const
Definition: AliEmcalJet.h:225
Bool_t TestFlavourTag(Int_t tag) const
Definition: AliEmcalJet.h:248
Double32_t fAreaEmc
Area on EMCAL surface (determined by ghosts in EMCal acceptance)
Definition: AliEmcalJet.h:284
UShort_t Nn() const
Definition: AliEmcalJet.h:142
AliVCluster * Cluster(Int_t idx) const
Double_t MaxClusterPt() const
Definition: AliEmcalJet.h:147
Represent a jet reconstructed using the EMCal jet framework.
Definition: AliEmcalJet.h:44
Int_t fNn
Number of neutral constituents.
Definition: AliEmcalJet.h:293
void SetNumberOfTracks(Int_t n)
Definition: AliEmcalJet.h:195
Bool_t PxPyPz(Double_t p[3]) const
Definition: AliEmcalJet.h:104
unsigned short UShort_t
Definition: External.C:28
Double_t Pz() const
Definition: AliEmcalJet.h:101
const std::vector< TLorentzVector > GetGhosts() const
Definition: AliEmcalJet.h:238
void SetTrigger(UInt_t trigger)
Definition: AliEmcalJet.h:213
void SetFlavour(Int_t flavour)
Definition: AliEmcalJet.h:247
const char Option_t
Definition: External.C:48
UInt_t fTriggers
! Triggers that the jet might have fired (AliVEvent::EOfflineTriggerTypes)
Definition: AliEmcalJet.h:308
Bool_t fAxisInEmcal
Whether the jet axis is inside the EMCAL acceptance.
Definition: AliEmcalJet.h:285
void SortConstituents()
Int_t fNEmc
Number of constituents in EMCAL acceptance.
Definition: AliEmcalJet.h:297
void CreateShapeProperties()
Definition: AliEmcalJet.h:258
Generic signal 1.
Definition: AliEmcalJet.h:80
Generic background 2.
Definition: AliEmcalJet.h:83
Double32_t fPhi
Jet axis azimuthal angle.
Definition: AliEmcalJet.h:270
DCal fiducial acceptance (each eta, phi edge narrowed by jet R)
Definition: AliEmcalJet.h:65
AliEmcalJet * fClosestJets[2]
! If this is MC it contains the two closest detector level jets in order of distance and viceversa ...
Definition: AliEmcalJet.h:300
bool Bool_t
Definition: External.C:53
AliEmcalJetShapeProperties * GetShapeProperties() const
Definition: AliEmcalJet.h:256
void AddTrigger(UInt_t trigger)
Definition: AliEmcalJet.h:214
DCal acceptance – spans ONLY DCal (no PHOS or gap)
Definition: AliEmcalJet.h:66
Double_t NEF() const
Definition: AliEmcalJet.h:141
void SetAreaPhi(Double_t a)
Definition: AliEmcalJet.h:187
void SetPtEmc(Double_t pt)
Definition: AliEmcalJet.h:200
Bool_t IsTriggerJet(UInt_t trigger=AliVEvent::kEMCEJE) const
Definition: AliEmcalJet.h:212
Double32_t fNEF
Jet Neutral Energy Fraction.
Definition: AliEmcalJet.h:274
Double32_t fClosestJetsDist[2]
! Distance from the two closest jets
Definition: AliEmcalJet.h:301
DCal fiducial acceptance (each eta, phi edge narrowed by jet R)
Definition: AliEmcalJet.h:67
UShort_t fMatchingType
! Matching type
Definition: AliEmcalJet.h:303
Double_t M() const
Definition: AliEmcalJet.h:113
EMCal fiducial acceptance (each eta, phi edge narrowed by jet R)
Definition: AliEmcalJet.h:63
Int_t Compare(const TObject *obj) const
Double_t MaxPartPt() const
Definition: AliEmcalJet.h:149
Generic background 3.
Definition: AliEmcalJet.h:84
void SetMatchedToClosest(UShort_t m)
Definition: AliEmcalJet.h:220
void SetNumberOfClusters(Int_t n)
Definition: AliEmcalJet.h:194
std::vector< int > GetPtSortedTrackConstituentIndexes(TClonesArray *tracks) const
Double_t Zv() const
Definition: AliEmcalJet.h:107
void SetAxisInEmcal(Bool_t b)
Definition: AliEmcalJet.h:190
TObjArray * fFlavourTracks
Definition: AliEmcalJet.h:287
void SetAreaEmc(Double_t a)
Definition: AliEmcalJet.h:189
void PrintConstituents(TClonesArray *tracks, TClonesArray *clusters) const
void SetNumberOfNeutrals(Int_t n)
Definition: AliEmcalJet.h:197