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