8 #include <TLorentzVector.h>
10 #include <TClonesArray.h>
13 #include "AliVParticle.h"
14 #include "AliVCluster.h"
15 #include "AliVEvent.h"
17 class AliEmcalJet : public AliVParticle
32 AliEmcalJet(Double_t px, Double_t py, Double_t pz);
33 AliEmcalJet(Double_t pt, Double_t eta, Double_t phi, Double_t m);
34 AliEmcalJet(const AliEmcalJet &jet);
35 AliEmcalJet& operator=(const AliEmcalJet &jet);
37 Double_t Px() const { return fPt*TMath::Cos(fPhi); }
38 Double_t Py() const { return fPt*TMath::Sin(fPhi); }
39 Double_t Pz() const { return fPt*TMath::SinH(fEta); }
40 Double_t Pt() const { return fPt; }
41 Double_t P() const { return fPt*TMath::CosH(fEta); }
42 Bool_t PxPyPz(Double_t p[3]) const { p[0]=Px();p[1]=Py();p[2]=Pz(); return 1; }
43 Double_t Xv() const { return 0.; }
44 Double_t Yv() const { return 0.; }
45 Double_t Zv() const { return 0.; }
46 Bool_t XvYvZv(Double_t x[3]) const { x[0]=0;x[1]=0;x[2]=0; return 1; }
47 Double_t OneOverPt() const { return 1./fPt; }
48 Double_t Phi() const { return fPhi; }
49 Double_t Theta() const { return 2*TMath::ATan(TMath::Exp(-fEta)); }
50 Double_t E() const { Double_t p=P(); return TMath::Sqrt(fM*fM+p*p); }
51 Double_t M() const { return fM; }
52 Double_t Eta() const { return fEta; }
53 Double_t Y() const { Double_t e = E(); Double_t pz = Pz(); return 0.5*TMath::Log((e+pz)/(e-pz)); }
54 Short_t Charge() const { return 0; }
55 Int_t GetLabel() const { return fLabel; }
56 Int_t PdgCode() const { return 0; }
57 const Double_t *PID() const { return 0; }
58 void GetMom(TLorentzVector &vec) const;
59 void Print(Option_t* option = "") const;
61 Double_t Area() const { return fArea; }
62 Double_t AreaPt() const { return fArea; }
63 Double_t AreaEta() const { return fAreaEta; }
64 Double_t AreaPhi() const { return fAreaPhi; }
65 Double_t AreaE() const { return fAreaE; }
66 Double_t AreaEmc() const { return fAreaEmc; }
67 Bool_t AxisInEmcal() const { return fAxisInEmcal; }
68 Int_t Compare(const TObject* obj) const;
69 Short_t ClusterAt(Int_t idx) const { return fClusterIDs.At(idx); }
70 AliVCluster *ClusterAt(Int_t idx, TClonesArray *ca) const { if (!ca) return 0; return dynamic_cast<AliVCluster*>(ca->At(ClusterAt(idx))); }
71 AliVCluster *GetLeadingCluster(TClonesArray *clusters) const;
72 UShort_t GetNumberOfClusters() const { return fClusterIDs.GetSize(); }
73 UShort_t GetNumberOfTracks() const { return fTrackIDs.GetSize(); }
74 UShort_t GetNumberOfConstituents() const { return GetNumberOfClusters()+GetNumberOfTracks(); }
75 Double_t FracEmcalArea() const { return fAreaEmc/fArea; }
76 Bool_t IsInsideEmcal() const { return (fAreaEmc/fArea>0.999); }
77 Bool_t IsInEmcal() const { return (fAreaEmc>0); }
78 Bool_t IsMC() const { return (Bool_t)(MCPt() > 0); }
79 Bool_t IsSortable() const { return kTRUE; }
80 Double_t MaxNeutralPt() const { return fMaxNPt; }
81 Double_t MaxChargedPt() const { return fMaxCPt; }
82 Double_t NEF() const { return fNEF; }
83 UShort_t Nn() const { return fNn; }
84 UShort_t Nch() const { return fNch; }
85 UShort_t N() const { return Nch()+Nn(); }
86 Int_t NEmc() const { return fNEmc; }
87 Double_t MCPt() const { return fMCPt; }
88 Double_t MaxClusterPt() const { return MaxNeutralPt(); }
89 Double_t MaxTrackPt() const { return MaxChargedPt(); }
90 Double_t MaxPartPt() const { return fMaxCPt < fMaxNPt ? fMaxNPt : fMaxCPt; }
91 Double_t PtEmc() const { return fPtEmc; }
92 Double_t PtSub() const { return fPtSub; }
93 Double_t PtSubVect() const { return fPtSubVect; }
94 Double_t PtSub(Double_t rho, Bool_t save = kFALSE);
95 Double_t PtSubVect(Double_t rho, Bool_t save = kFALSE);
96 TLorentzVector SubtractRhoVect(Double_t rho, Bool_t save = kFALSE);
97 Short_t TrackAt(Int_t idx) const { return fTrackIDs.At(idx); }
98 AliVParticle *TrackAt(Int_t idx, TClonesArray *ta) const { if (!ta) return 0; return dynamic_cast<AliVParticle*>(ta->At(TrackAt(idx))); }
99 AliVParticle *GetLeadingTrack(TClonesArray *tracks) const;
100 Int_t GetFlavour() const { return fFlavourTagging; }
102 void AddClusterAt(Int_t clus, Int_t idx){ fClusterIDs.AddAt(clus, idx); }
103 void AddFlavourTag(Int_t tag) { fFlavourTagging |= tag; }
104 void AddTrackAt(Int_t track, Int_t idx) { fTrackIDs.AddAt(track, idx); }
105 void Clear(Option_t */*option*/="") { fClusterIDs.Set(0); fTrackIDs.Set(0); fClosestJets[0] = 0; fClosestJets[1] = 0;
106 fClosestJetsDist[0] = 0; fClosestJetsDist[1] = 0; fMatched = 0; fPtSub = 0; }
107 Double_t DeltaR(const AliVParticle* part) const;
108 Double_t GetZ ( const Double_t trkPx, const Double_t trkPy, const Double_t trkPz ) const; // Get Z of constituent trk
109 Double_t GetZ ( const AliVParticle* trk ) const; // Get Z of constituent trk
110 Double_t GetXi ( const AliVParticle* trk ) const { return TMath::Log ( 1/GetZ (trk) ); } // Get Xi of constituent trk
111 Double_t GetXi ( const Double_t trkPx, const Double_t trkPy, const Double_t trkPz ) const { return TMath::Log ( 1/GetZ (trkPx, trkPy, trkPz ) ); } // Get Xi of constituent trk
113 void SetLabel(Int_t l) { fLabel = l; }
114 void SetArea(Double_t a) { fArea = a; }
115 void SetAreaEta(Double_t a) { fAreaEta = a; }
116 void SetAreaPhi(Double_t a) { fAreaPhi = TVector2::Phi_0_2pi(a); }
117 void SetAreaE(Double_t a) { fAreaE = a; }
118 void SetAreaEmc(Double_t a) { fAreaEmc = a; }
119 void SetAxisInEmcal(Bool_t b) { fAxisInEmcal = b; }
120 void SetFlavour(Int_t flavour) { fFlavourTagging = flavour; }
121 void SetMaxNeutralPt(Double32_t t) { fMaxNPt = t; }
122 void SetMaxChargedPt(Double32_t t) { fMaxCPt = t; }
123 void SetNEF(Double_t nef) { fNEF = nef; }
124 void SetNumberOfClusters(Int_t n) { fClusterIDs.Set(n); }
125 void SetNumberOfTracks(Int_t n) { fTrackIDs.Set(n); }
126 void SetNumberOfCharged(Int_t n) { fNch = n; }
127 void SetNumberOfNeutrals(Int_t n) { fNn = n; }
128 void SetMCPt(Double_t p) { fMCPt = p; }
129 void SortConstituents();
130 std::vector<int> SortConstituentsPt(TClonesArray *tracks) const;
131 void SetNEmc(Int_t n) { fNEmc = n; }
132 void SetPtEmc(Double_t pt) { fPtEmc = pt; }
133 void SetPtSub(Double_t ps) { fPtSub = ps; }
134 void SetPtSubVect(Double_t ps) { fPtSubVect = ps; }
135 Bool_t TestFlavourTag(Int_t tag) const { return (Bool_t)((tag & fFlavourTagging) !=0); }
138 Bool_t IsTriggerJet(UInt_t trigger=AliVEvent::kEMCEJE) const { return (Bool_t)((fTriggers & trigger) != 0); }
139 void SetTrigger(UInt_t trigger) { fTriggers = trigger; }
140 void AddTrigger(UInt_t trigger) { fTriggers |= trigger; }
143 void SetClosestJet(AliEmcalJet *j, Double_t d) { fClosestJets[0] = j; fClosestJetsDist[0] = d ; }
144 void SetSecondClosestJet(AliEmcalJet *j, Double_t d) { fClosestJets[1] = j; fClosestJetsDist[1] = d ; }
145 void SetMatchedToClosest(UShort_t m) { fMatched = 0; fMatchingType = m ; }
146 void SetMatchedToSecondClosest(UShort_t m) { fMatched = 1; fMatchingType = m ; }
147 void ResetMatching();
148 AliEmcalJet* ClosestJet() const { return fClosestJets[0] ; }
149 Double_t ClosestJetDistance() const { return fClosestJetsDist[0] ; }
150 AliEmcalJet* SecondClosestJet() const { return fClosestJets[1] ; }
151 Double_t SecondClosestJetDistance() const { return fClosestJetsDist[1] ; }
152 AliEmcalJet* MatchedJet() const { return fMatched < 2 ? fClosestJets[fMatched] : 0; }
153 UShort_t GetMatchingType() const { return fMatchingType ; }
155 void SetTaggedJet(AliEmcalJet *j) { fTaggedJet = j ; }
156 void SetTagStatus(Int_t i) { fTagStatus = i ; }
157 AliEmcalJet* GetTaggedJet() const { return fTaggedJet ; }
158 Int_t GetTagStatus() const { return fTagStatus ; }
160 //jet shape derivatives
162 void SetFirstDerivative(Double_t d) { fJetShapeMassFirstDer = d ; }
163 void SetSecondDerivative(Double_t d) { fJetShapeMassSecondDer = d ; }
164 void SetFirstOrderSubtracted(Double_t d) { fJetShapeMassFirstSub = d ; }
165 void SetSecondOrderSubtracted(Double_t d) { fJetShapeMassSecondSub = d ; }
166 Double_t GetFirstDerivative() const { return fJetShapeMassFirstDer ; }
167 Double_t GetSecondDerivative() const { return fJetShapeMassSecondDer ; }
168 Double_t GetFirstOrderSubtracted() const { return fJetShapeMassFirstSub ; }
169 Double_t GetSecondOrderSubtracted() const { return fJetShapeMassSecondSub ; }
171 //jet structure function
172 TArrayF GetGRNumerator() const { return fGRNumerator ; }
173 TArrayF GetGRDenominator() const { return fGRDenominator ; }
174 TArrayF GetGRNumeratorSub() const { return fGRNumeratorSub ; }
175 TArrayF GetGRDenominatorSub() const { return fGRDenominatorSub ; }
176 void AddGRNumAt(Float_t num, Int_t idx) { fGRNumerator.AddAt(num, idx) ; }
177 void AddGRDenAt(Float_t den, Int_t idx) { fGRDenominator.AddAt(den, idx) ; }
178 void SetGRNumSize(UInt_t s) { fGRNumerator.Set(s) ; }
179 void SetGRDenSize(UInt_t s) { fGRDenominator.Set(s) ; }
181 void AddGRNumSubAt(Float_t num, Int_t idx) { fGRNumeratorSub.AddAt(num, idx) ; }
182 void AddGRDenSubAt(Float_t den, Int_t idx) { fGRDenominatorSub.AddAt(den, idx) ; }
183 void SetGRNumSubSize(UInt_t s) { fGRNumeratorSub.Set(s) ; }
184 void SetGRDenSubSize(UInt_t s) { fGRDenominatorSub.Set(s) ; }
188 void SetFirstDerivativeAngularity(Double_t d) { fJetShapeAngularityFirstDer = d ; }
189 void SetSecondDerivativeAngularity(Double_t d) { fJetShapeAngularitySecondDer = d ; }
190 void SetFirstOrderSubtractedAngularity(Double_t d) { fJetShapeAngularityFirstSub = d ; }
191 void SetSecondOrderSubtractedAngularity(Double_t d) { fJetShapeAngularitySecondSub = d ; }
192 Double_t GetFirstDerivativeAngularity() const { return fJetShapeAngularityFirstDer ; }
193 Double_t GetSecondDerivativeAngularity() const { return fJetShapeAngularitySecondDer ; }
194 Double_t GetFirstOrderSubtractedAngularity() const { return fJetShapeAngularityFirstSub ; }
195 Double_t GetSecondOrderSubtractedAngularity() const { return fJetShapeAngularitySecondSub ; }
198 void SetFirstDerivativepTD(Double_t d) { fJetShapepTDFirstDer = d ; }
199 void SetSecondDerivativepTD(Double_t d) { fJetShapepTDSecondDer = d ; }
200 void SetFirstOrderSubtractedpTD(Double_t d) { fJetShapepTDFirstSub = d ; }
201 void SetSecondOrderSubtractedpTD(Double_t d) { fJetShapepTDSecondSub = d ; }
202 Double_t GetFirstDerivativepTD() const { return fJetShapepTDFirstDer ; }
203 Double_t GetSecondDerivativepTD() const { return fJetShapepTDSecondDer ; }
204 Double_t GetFirstOrderSubtractedpTD() const { return fJetShapepTDFirstSub ; }
205 Double_t GetSecondOrderSubtractedpTD() const { return fJetShapepTDSecondSub ; }
208 void SetFirstDerivativeCircularity(Double_t d) { fJetShapeCircularityFirstDer = d ; }
209 void SetSecondDerivativeCircularity(Double_t d) { fJetShapeCircularitySecondDer = d ; }
210 void SetFirstOrderSubtractedCircularity(Double_t d) { fJetShapeCircularityFirstSub = d ; }
211 void SetSecondOrderSubtractedCircularity(Double_t d) { fJetShapeCircularitySecondSub = d ; }
212 Double_t GetFirstDerivativeCircularity() const { return fJetShapeCircularityFirstDer ; }
213 Double_t GetSecondDerivativeCircularity() const { return fJetShapeCircularitySecondDer ; }
214 Double_t GetFirstOrderSubtractedCircularity() const { return fJetShapeCircularityFirstSub ; }
215 Double_t GetSecondOrderSubtractedCircularity() const { return fJetShapeCircularitySecondSub ; }
218 void SetFirstDerivativeSigma2(Double_t d) { fJetShapeSigma2FirstDer = d ; }
219 void SetSecondDerivativeSigma2(Double_t d) { fJetShapeSigma2SecondDer = d ; }
220 void SetFirstOrderSubtractedSigma2(Double_t d) { fJetShapeSigma2FirstSub = d ; }
221 void SetSecondOrderSubtractedSigma2(Double_t d) { fJetShapeSigma2SecondSub = d ; }
222 Double_t GetFirstDerivativeSigma2() const { return fJetShapeSigma2FirstDer ; }
223 Double_t GetSecondDerivativeSigma2() const { return fJetShapeSigma2SecondDer ; }
224 Double_t GetFirstOrderSubtractedSigma2() const { return fJetShapeSigma2FirstSub ; }
225 Double_t GetSecondOrderSubtractedSigma2() const { return fJetShapeSigma2SecondSub ; }
228 //number of contituents
229 void SetFirstDerivativeConstituent(Double_t d) { fJetShapeConstituentFirstDer = d ; }
230 void SetSecondDerivativeConstituent(Double_t d) { fJetShapeConstituentSecondDer = d ; }
231 void SetFirstOrderSubtractedConstituent(Double_t d) { fJetShapeConstituentFirstSub = d ; }
232 void SetSecondOrderSubtractedConstituent(Double_t d) { fJetShapeConstituentSecondSub = d ; }
233 Double_t GetFirstDerivativeConstituent() const { return fJetShapeConstituentFirstDer ; }
234 Double_t GetSecondDerivativeConstituent() const { return fJetShapeConstituentSecondDer ; }
235 Double_t GetFirstOrderSubtractedConstituent() const { return fJetShapeConstituentFirstSub ; }
236 Double_t GetSecondOrderSubtractedConstituent() const { return fJetShapeConstituentSecondSub ; }
238 //leading minus subleading constituent
239 void SetFirstDerivativeLeSub(Double_t d) { fJetShapeLeSubFirstDer = d ; }
240 void SetSecondDerivativeLeSub(Double_t d) { fJetShapeLeSubSecondDer = d ; }
241 void SetFirstOrderSubtractedLeSub(Double_t d) { fJetShapeLeSubFirstSub = d ; }
242 void SetSecondOrderSubtractedLeSub(Double_t d) { fJetShapeLeSubSecondSub = d ; }
243 Double_t GetFirstDerivativeLeSub() const { return fJetShapeLeSubFirstDer ; }
244 Double_t GetSecondDerivativeLeSub() const { return fJetShapeLeSubSecondDer ; }
245 Double_t GetFirstOrderSubtractedLeSub() const { return fJetShapeLeSubFirstSub ; }
246 Double_t GetSecondOrderSubtractedLeSub() const { return fJetShapeLeSubSecondSub ; }
249 Double32_t fPt; //[0,0,12] pt
250 Double32_t fEta; //[-1,1,12] eta
251 Double32_t fPhi; //[0,6.3,12] phi
252 Double32_t fM; //[0,0,8] mass
253 Double32_t fNEF; //[0,1,8] neutral energy fraction
254 Double32_t fArea; //[0,0,12] area (transverse)
255 Double32_t fAreaEta; //[0,0,12] area eta
256 Double32_t fAreaPhi; //[0,0,12] area phi
257 Double32_t fAreaE; //[0,0,12] temporal area component
258 Double32_t fAreaEmc; //[0,0,12] area on EMCAL surface (determined from ghosts)
259 Bool_t fAxisInEmcal; // =true if jet axis inside EMCAL acceptance
260 Int_t fFlavourTagging; // tag jet with a falvour, bit 0 = no tag; bit 1= Dstar; bit 2 = D0
261 Double32_t fMaxCPt; //[0,0,12] pt of maximum charged constituent
262 Double32_t fMaxNPt; //[0,0,12] pt of maximum neutral constituent
263 Double32_t fMCPt; // pt from MC particles contributing to the jet
264 Int_t fNn; // number of neutral constituents
265 Int_t fNch; // number of charged constituents
266 Double32_t fPtEmc; //[0,0,12] pt in EMCAL acceptance
267 Int_t fNEmc; // number of constituents in EMCAL acceptance
268 TArrayI fClusterIDs; // array containing ids of cluster constituents
269 TArrayI fTrackIDs; // array containing ids of track constituents
270 AliEmcalJet *fClosestJets[2]; //! if this is MC it contains the two closest detector level jets in order of distance and viceversa
271 Double32_t fClosestJetsDist[2]; //! distance to closest jets (see above)
272 UShort_t fMatched; //! 0,1 if it is matched with one of the closest jets; 2 if it is not matched
273 UShort_t fMatchingType; //! matching type
274 AliEmcalJet *fTaggedJet; //! jet tagged to this jet
275 Int_t fTagStatus; //! status of tagging -1: NA 0: not tagged 1: tagged
276 Double_t fPtSub; //! background subtracted pt (not stored set from outside)
277 Double_t fPtSubVect; //! background vector subtracted pt (not stored set from outside)
278 UInt_t fTriggers; //! triggers that the jet might have fired (AliVEvent::EOfflineTriggerTypes)
280 Double_t fJetShapeMassFirstDer; //! result from shape derivatives for jet mass: 1st derivative
281 Double_t fJetShapeMassSecondDer; //! result from shape derivatives for jet mass: 2nd derivative
282 Double_t fJetShapeMassFirstSub; //! result from shape derivatives for jet mass: 1st order subtracted
283 Double_t fJetShapeMassSecondSub; //! result from shape derivatives for jet mass: 2nd order subtracted
284 Int_t fLabel; // label to inclusive jet for constituent subtracted jet
286 TArrayF fGRNumerator; //! array with angular structure function numerator
287 TArrayF fGRDenominator; //! array with angular structure function denominator
288 TArrayF fGRNumeratorSub; //! array with angular structure function numerator
289 TArrayF fGRDenominatorSub; //! array with angular structure function denominator
291 Double_t fJetShapeAngularityFirstDer; //! result from shape derivatives for jet Angularity: 1st derivative
292 Double_t fJetShapeAngularitySecondDer; //! result from shape derivatives for jet Angularity: 2nd derivative
293 Double_t fJetShapeAngularityFirstSub; //! result from shape derivatives for jet Angularity: 1st order subtracted
294 Double_t fJetShapeAngularitySecondSub; //! result from shape derivatives for jet Angularity: 2nd order subtracted
296 Double_t fJetShapepTDFirstDer; //! result from shape derivatives for jet pTD: 1st derivative
297 Double_t fJetShapepTDSecondDer; //! result from shape derivatives for jet pTD: 2nd derivative
298 Double_t fJetShapepTDFirstSub; //! result from shape derivatives for jet pTD: 1st order subtracted
299 Double_t fJetShapepTDSecondSub; //! result from shape derivatives for jet pTD: 2nd order subtracted
301 Double_t fJetShapeCircularityFirstDer; //! result from shape derivatives for jet circularity: 1st derivative
302 Double_t fJetShapeCircularitySecondDer; //! result from shape derivatives for jet circularity: 2nd derivative
303 Double_t fJetShapeCircularityFirstSub; //! result from shape derivatives for jet circularity: 1st order subtracted
304 Double_t fJetShapeCircularitySecondSub; //! result from shape derivatives for jetcircularity: 2nd order subtracted
306 Double_t fJetShapeSigma2FirstDer; //! result from shape derivatives for jet sigma2: 1st derivative
307 Double_t fJetShapeSigma2SecondDer; //! result from shape derivatives for jet sigma2: 2nd derivative
308 Double_t fJetShapeSigma2FirstSub; //! result from shape derivatives for jet sigma2: 1st order subtracted
309 Double_t fJetShapeSigma2SecondSub; //! result from shape derivatives for jetsigma2: 2nd order subtracted
311 Double_t fJetShapeConstituentFirstDer; //! result from shape derivatives for jet const: 1st derivative
312 Double_t fJetShapeConstituentSecondDer; //! result from shape derivatives for jet const: 2nd derivative
313 Double_t fJetShapeConstituentFirstSub; //! result from shape derivatives for jet const: 1st order subtracted
314 Double_t fJetShapeConstituentSecondSub; //! result from shape derivatives for jet const: 2nd order subtracted
316 Double_t fJetShapeLeSubFirstDer; //! result from shape derivatives for jet LeSub: 1st derivative
317 Double_t fJetShapeLeSubSecondDer; //! result from shape derivatives for jet LeSub: 2nd derivative
318 Double_t fJetShapeLeSubFirstSub; //! result from shape derivatives for jet LeSub: 1st order subtracted
319 Double_t fJetShapeLeSubSecondSub; //! result from shape derivatives for jet LeSub: 2nd order subtracted
323 { // sort in decreasing order
324 // first value of the pair is Pt and the second is entry index
325 bool operator () (const std::pair<Double_t, Int_t>& p1, const std::pair<Double_t, Int_t>& p2) { return p1.first > p2.first ; }
328 ClassDef(AliEmcalJet,16) // Emcal jet class in cylindrical coordinates