]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/AliEmcalJet.h
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / AliEmcalJet.h
1 #ifndef AliEmcalJet_H
2 #define AliEmcalJet_H
3
4 #include <vector>
5 #include <algorithm>
6 #include <utility>
7 #include <TArrayS.h>
8 #include <TLorentzVector.h>
9 #include <TMath.h>
10 #include <TClonesArray.h>
11 #include <TVector2.h>
12
13 #include "AliVParticle.h"
14 #include "AliVCluster.h"
15 #include "AliVEvent.h"
16
17 class AliEmcalJet : public AliVParticle
18 {
19  public:
20      enum EFlavourTag{
21        kDStar = 1<<0,
22        kD0 = 1<<1,
23        kSig1 = 1<<2,
24        kSig2 = 1<<3,
25        kBckgrd1 = 1<<4,
26        kBckgrd2 = 1<<5,
27        kBckgrd3 = 1<<6
28        //.....
29     }; 
30  
31   AliEmcalJet();
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);
36
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(M()*M()+p*p); }
51   Double_t          M()                          const { return fM; }
52   Double_t          Eta()                        const { return fEta;    }
53   Double_t          Y()                          const { 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;
60
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          AreaEmc()                    const { return fAreaEmc;                  }
66   Bool_t            AxisInEmcal()                const { return fAxisInEmcal;              }
67   Int_t             Compare(const TObject* obj)  const;
68   Short_t           ClusterAt(Int_t idx)         const { return fClusterIDs.At(idx);       }
69   AliVCluster      *ClusterAt(Int_t idx, TClonesArray *ca)  const { if (!ca) return 0; return dynamic_cast<AliVCluster*>(ca->At(ClusterAt(idx))); }
70   AliVCluster      *GetLeadingCluster(TClonesArray *clusters) const;
71   UShort_t          GetNumberOfClusters()        const { return fClusterIDs.GetSize();     }
72   UShort_t          GetNumberOfTracks()          const { return fTrackIDs.GetSize();       }
73   UShort_t          GetNumberOfConstituents()    const { return GetNumberOfClusters()+GetNumberOfTracks();       }
74   Double_t          FracEmcalArea()              const { return fAreaEmc/fArea;            }
75   Bool_t            IsInsideEmcal()              const { return (fAreaEmc/fArea>0.999);    }
76   Bool_t            IsInEmcal()                  const { return (fAreaEmc>0);              }
77   Bool_t            IsMC()                       const { return (Bool_t)(MCPt() > 0);      }
78   Bool_t            IsSortable()                 const { return kTRUE;                     }
79   Double_t          MaxNeutralPt()               const { return fMaxNPt;                   }
80   Double_t          MaxChargedPt()               const { return fMaxCPt;                   }
81   Double_t          NEF()                        const { return fNEF;                      }
82   UShort_t          Nn()                         const { return fNn;                       }
83   UShort_t          Nch()                        const { return fNch;                      }
84   UShort_t          N()                          const { return Nch()+Nn();                }
85   Int_t             NEmc()                       const { return fNEmc;                     }
86   Double_t          MCPt()                       const { return fMCPt;                     }
87   Double_t          MaxClusterPt()               const { return MaxNeutralPt();            }
88   Double_t          MaxTrackPt()                 const { return MaxChargedPt();            }
89   Double_t          MaxPartPt()                  const { return fMaxCPt < fMaxNPt ? fMaxNPt : fMaxCPt;     }
90   Double_t          PtEmc()                      const { return fPtEmc;                    }
91   Double_t          PtSub()                      const { return fPtSub;                    }
92   Double_t          PtSub(Double_t rho)          const { return fPt - fArea*rho;           }
93   Double_t          PtSubVect(Double_t rho)      const;
94   Short_t           TrackAt(Int_t idx)           const { return fTrackIDs.At(idx);         }
95   AliVParticle     *TrackAt(Int_t idx, TClonesArray *ta)  const { if (!ta) return 0; return dynamic_cast<AliVParticle*>(ta->At(TrackAt(idx))); } 
96   AliVParticle     *GetLeadingTrack(TClonesArray *tracks) const;
97   Int_t             GetFlavour()                 const { return fFlavourTagging;           } 
98   
99   void              AddClusterAt(Int_t clus, Int_t idx){ fClusterIDs.AddAt(clus, idx);     }
100   void              AddFlavourTag(Int_t tag)           { fFlavourTagging |= tag; }
101   void              AddTrackAt(Int_t track, Int_t idx) { fTrackIDs.AddAt(track, idx);      }
102   void              Clear(Option_t */*option*/="")     { fClusterIDs.Set(0); fTrackIDs.Set(0); fClosestJets[0] = 0; fClosestJets[1] = 0; 
103                                                          fClosestJetsDist[0] = 0; fClosestJetsDist[1] = 0; fMatched = 0; fPtSub = 0; }
104   Double_t          DeltaR(const AliVParticle* part) const;
105   Double_t          GetZ ( const Double_t trkPx, const Double_t trkPy, const Double_t trkPz ) const; // Get Z of constituent trk
106   Double_t          GetZ ( const AliVParticle* trk )       const; // Get Z of constituent trk
107   Double_t          GetXi ( const AliVParticle* trk )      const { return TMath::Log ( 1/GetZ (trk) ); } // Get Xi of constituent trk
108   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
109
110   void              SetLabel(Int_t l)                  { fLabel = l;                       }
111   void              SetArea(Double_t a)                { fArea    = a;                     }
112   void              SetAreaEta(Double_t a)             { fAreaEta = a;                     }
113   void              SetAreaPhi(Double_t a)             { fAreaPhi = a;                     }
114   void              SetAreaEmc(Double_t a)             { fAreaEmc = a;                     }
115   void              SetAxisInEmcal(Bool_t b)           { fAxisInEmcal = b;                 }
116   void              SetFlavour(Int_t flavour)          { fFlavourTagging = flavour;        }
117   void              SetMaxNeutralPt(Double32_t t)      { fMaxNPt  = t;                     }
118   void              SetMaxChargedPt(Double32_t t)      { fMaxCPt  = t;                     }
119   void              SetNEF(Double_t nef)               { fNEF     = nef;                   }
120   void              SetNumberOfClusters(Int_t n)       { fClusterIDs.Set(n);               }
121   void              SetNumberOfTracks(Int_t n)         { fTrackIDs.Set(n);                 }
122   void              SetNumberOfCharged(Int_t n)        { fNch = n;                         }
123   void              SetNumberOfNeutrals(Int_t n)       { fNn = n;                          }
124   void              SetMCPt(Double_t p)                { fMCPt = p;                        }
125   void              SortConstituents();
126   std::vector<int>  SortConstituentsPt(TClonesArray *tracks) const;
127   void              SetNEmc(Int_t n)                   { fNEmc           = n;              }
128   void              SetPtEmc(Double_t pt)              { fPtEmc          = pt;             }
129   void              SetPtSub(Double_t ps)              { fPtSub          = ps;             } 
130   void              SetPtSubVect(Double_t ps)          { fPtVectSub      = ps;             }
131   Bool_t            TestFlavourTag(Int_t tag)          { return (Bool_t)((tag & fFlavourTagging) !=0); }
132
133   // Trigger
134   Bool_t            IsTriggerJet(UInt_t trigger=AliVEvent::kEMCEJE) const   { return (Bool_t)((fTriggers & trigger) != 0); }
135   void              SetTrigger(UInt_t trigger)                              { fTriggers  = trigger;                        }
136   void              AddTrigger(UInt_t trigger)                              { fTriggers |= trigger;                        }
137
138   // Matching
139   void              SetClosestJet(AliEmcalJet *j, Double_t d)       { fClosestJets[0] = j; fClosestJetsDist[0] = d    ; }
140   void              SetSecondClosestJet(AliEmcalJet *j, Double_t d) { fClosestJets[1] = j; fClosestJetsDist[1] = d    ; }
141   void              SetMatchedToClosest(UShort_t m)                 { fMatched        = 0; fMatchingType       = m    ; }
142   void              SetMatchedToSecondClosest(UShort_t m)           { fMatched        = 1; fMatchingType       = m    ; }
143   void              ResetMatching();
144   AliEmcalJet*      ClosestJet()                              const { return fClosestJets[0]                          ; }
145   Double_t          ClosestJetDistance()                      const { return fClosestJetsDist[0]                      ; }
146   AliEmcalJet*      SecondClosestJet()                        const { return fClosestJets[1]                          ; }
147   Double_t          SecondClosestJetDistance()                const { return fClosestJetsDist[1]                      ; }
148   AliEmcalJet*      MatchedJet()                              const { return fMatched < 2 ? fClosestJets[fMatched] : 0; }
149   UShort_t          GetMatchingType()                         const { return fMatchingType                            ; }
150
151   void              SetTaggedJet(AliEmcalJet *j)                    { fTaggedJet = j                                  ; }
152   void              SetTagStatus(Int_t i)                           { fTagStatus = i                                  ; }
153   AliEmcalJet*      GetTaggedJet()                            const { return fTaggedJet                               ; }
154   Int_t             GetTagStatus()                            const { return fTagStatus                               ; }
155
156   //jet shape derivatives
157   //jet mass
158   void              SetFirstDerivative(Double_t d)                  { fJetShapeMassFirstDer = d                       ; }
159   void              SetSecondDerivative(Double_t d)                 { fJetShapeMassSecondDer = d                      ; }
160   void              SetFirstOrderSubtracted(Double_t d)             { fJetShapeMassFirstSub = d                       ; }
161   void              SetSecondOrderSubtracted(Double_t d)            { fJetShapeMassSecondSub = d                      ; }
162   Double_t          GetFirstDerivative()                      const { return fJetShapeMassFirstDer                    ; }
163   Double_t          GetSecondDerivative()                     const { return fJetShapeMassSecondDer                   ; }
164   Double_t          GetFirstOrderSubtracted()                 const { return fJetShapeMassFirstSub                    ; }
165   Double_t          GetSecondOrderSubtracted()                const { return fJetShapeMassSecondSub                   ; }
166   
167   //jet structure function
168   TArrayF           GetGRNumerator()                          const { return fGRNumerator                             ; }
169   TArrayF           GetGRDenominator()                        const { return fGRDenominator                           ; }
170   TArrayF           GetGRNumeratorSub()                       const { return fGRNumeratorSub                          ; }
171   TArrayF           GetGRDenominatorSub()                     const { return fGRDenominatorSub                        ; }
172   void              AddGRNumAt(Float_t num, Int_t idx)              { fGRNumerator.AddAt(num, idx)                    ; }
173   void              AddGRDenAt(Float_t den, Int_t idx)              { fGRDenominator.AddAt(den, idx)                  ; }
174   void              SetGRNumSize(UInt_t s)                          { fGRNumerator.Set(s)                             ; }
175   void              SetGRDenSize(UInt_t s)                          { fGRDenominator.Set(s)                           ; }
176
177   void              AddGRNumSubAt(Float_t num, Int_t idx)           { fGRNumeratorSub.AddAt(num, idx)                 ; }
178   void              AddGRDenSubAt(Float_t den, Int_t idx)           { fGRDenominatorSub.AddAt(den, idx)               ; }
179   void              SetGRNumSubSize(UInt_t s)                       { fGRNumeratorSub.Set(s)                          ; }
180   void              SetGRDenSubSize(UInt_t s)                       { fGRDenominatorSub.Set(s)                        ; }
181   void              PrintGR();
182
183   //Angularity
184   void              SetFirstDerivativeAngularity(Double_t d)        { fJetShapeAngularityFirstDer = d                 ; }
185   void              SetSecondDerivativeAngularity(Double_t d)       { fJetShapeAngularitySecondDer = d                ; }
186   void              SetFirstOrderSubtractedAngularity(Double_t d)   { fJetShapeAngularityFirstSub = d                 ; }
187   void              SetSecondOrderSubtractedAngularity(Double_t d)  { fJetShapeAngularitySecondSub = d                ; }
188   Double_t          GetFirstDerivativeAngularity()            const { return fJetShapeAngularityFirstDer              ; }
189   Double_t          GetSecondDerivativeAngularity()           const { return fJetShapeAngularitySecondDer             ; }
190   Double_t          GetFirstOrderSubtractedAngularity()       const { return fJetShapeAngularityFirstSub              ; }
191   Double_t          GetSecondOrderSubtractedAngularity()      const { return fJetShapeAngularitySecondSub             ; }
192
193   //pTD
194   void              SetFirstDerivativepTD(Double_t d)               { fJetShapepTDFirstDer = d                        ; }
195   void              SetSecondDerivativepTD(Double_t d)              { fJetShapepTDSecondDer = d                       ; }
196   void              SetFirstOrderSubtractedpTD(Double_t d)          { fJetShapepTDFirstSub = d                        ; }
197   void              SetSecondOrderSubtractedpTD(Double_t d)         { fJetShapepTDSecondSub = d                       ; }
198   Double_t          GetFirstDerivativepTD()                   const { return fJetShapepTDFirstDer                     ; }
199   Double_t          GetSecondDerivativepTD()                  const { return fJetShapepTDSecondDer                    ; }
200   Double_t          GetFirstOrderSubtractedpTD()              const { return fJetShapepTDFirstSub                     ; }
201   Double_t          GetSecondOrderSubtractedpTD()             const { return fJetShapepTDSecondSub                    ; }
202  
203   //Circularity
204   void              SetFirstDerivativeCircularity(Double_t d)       { fJetShapeCircularityFirstDer = d                ; }
205   void              SetSecondDerivativeCircularity(Double_t d)      { fJetShapeCircularitySecondDer = d               ; }
206   void              SetFirstOrderSubtractedCircularity(Double_t d)  { fJetShapeCircularityFirstSub = d                ; }
207   void              SetSecondOrderSubtractedCircularity(Double_t d) { fJetShapeCircularitySecondSub = d               ; }
208   Double_t          GetFirstDerivativeCircularity()           const { return fJetShapeCircularityFirstDer             ; }
209   Double_t          GetSecondDerivativeCircularity()          const { return fJetShapeCircularitySecondDer            ; }
210   Double_t          GetFirstOrderSubtractedCircularity()      const { return fJetShapeCircularityFirstSub             ; }
211   Double_t          GetSecondOrderSubtractedCircularity()     const { return fJetShapeCircularitySecondSub            ; }
212
213   //number of contituents
214   void              SetFirstDerivativeConstituent(Double_t d)       { fJetShapeConstituentFirstDer = d                ; }
215   void              SetSecondDerivativeConstituent(Double_t d)      { fJetShapeConstituentSecondDer = d               ; }
216   void              SetFirstOrderSubtractedConstituent(Double_t d)  { fJetShapeConstituentFirstSub = d                ; }
217   void              SetSecondOrderSubtractedConstituent(Double_t d) { fJetShapeConstituentSecondSub = d               ; }
218   Double_t          GetFirstDerivativeConstituent()           const { return fJetShapeConstituentFirstDer             ; }
219   Double_t          GetSecondDerivativeConstituent()          const { return fJetShapeConstituentSecondDer            ; }
220   Double_t          GetFirstOrderSubtractedConstituent()      const { return fJetShapeConstituentFirstSub             ; }
221   Double_t          GetSecondOrderSubtractedConstituent()     const { return fJetShapeConstituentSecondSub            ; }
222
223   //leading minus subleading constituent
224   void              SetFirstDerivativeLeSub(Double_t d)             { fJetShapeLeSubFirstDer = d                      ; }
225   void              SetSecondDerivativeLeSub(Double_t d)            { fJetShapeLeSubSecondDer = d                     ; }
226   void              SetFirstOrderSubtractedLeSub(Double_t d)        { fJetShapeLeSubFirstSub = d                      ; }
227   void              SetSecondOrderSubtractedLeSub(Double_t d)       { fJetShapeLeSubSecondSub = d                     ; }
228   Double_t          GetFirstDerivativeLeSub()                 const { return fJetShapeLeSubFirstDer                   ; }
229   Double_t          GetSecondDerivativeLeSub()                const { return fJetShapeLeSubSecondDer                  ; }
230   Double_t          GetFirstOrderSubtractedLeSub()            const { return fJetShapeLeSubFirstSub                   ; }
231   Double_t          GetSecondOrderSubtractedLeSub()           const { return fJetShapeLeSubSecondSub                  ; }
232   
233  protected:
234   Double32_t        fPt;                  //[0,0,12]   pt 
235   Double32_t        fEta;                 //[-1,1,12]  eta
236   Double32_t        fPhi;                 //[0,6.3,12] phi
237   Double32_t        fM;                   //[0,0,8]    mass
238   Double32_t        fNEF;                 //[0,1,8]    neutral energy fraction
239   Double32_t        fArea;                //[0,0,12]   area
240   Double32_t        fAreaEta;             //[0,0,12]   area eta
241   Double32_t        fAreaPhi;             //[0,0,12]   area phi
242   Double32_t        fAreaEmc;             //[0,0,12]   area on EMCAL surface (determined from ghosts)
243   Bool_t            fAxisInEmcal;         //           =true if jet axis inside EMCAL acceptance
244   Int_t             fFlavourTagging;      // tag jet with a falvour, bit 0 = no tag; bit 1= Dstar; bit 2 = D0
245   Double32_t        fMaxCPt;              //[0,0,12]   pt of maximum charged constituent
246   Double32_t        fMaxNPt;              //[0,0,12]   pt of maximum neutral constituent
247   Double32_t        fMCPt;                //           pt from MC particles contributing to the jet
248   Int_t             fNn;                  //           number of neutral constituents
249   Int_t             fNch;                 //           number of charged constituents
250   Double32_t        fPtEmc;               //[0,0,12]   pt in EMCAL acceptance
251   Int_t             fNEmc;                //           number of constituents in EMCAL acceptance
252   TArrayI           fClusterIDs;          //           array containing ids of cluster constituents  
253   TArrayI           fTrackIDs;            //           array containing ids of track constituents   
254   AliEmcalJet      *fClosestJets[2];      //!          if this is MC it contains the two closest detector level jets in order of distance and viceversa
255   Double32_t        fClosestJetsDist[2];  //!          distance to closest jets (see above)
256   UShort_t          fMatched;             //!          0,1 if it is matched with one of the closest jets; 2 if it is not matched
257   UShort_t          fMatchingType;        //!          matching type
258   AliEmcalJet      *fTaggedJet;           //!          jet tagged to this jet
259   Int_t             fTagStatus;           //!          status of tagging -1: NA 0: not tagged 1: tagged
260   Double_t          fPtSub;               //!          background subtracted pt (not stored set from outside) 
261   Double_t          fPtVectSub;           //!          background vector subtracted pt (not stored set from outside) 
262   UInt_t            fTriggers;            //!          triggers that the jet might have fired (AliVEvent::EOfflineTriggerTypes)
263
264   Double_t          fJetShapeMassFirstDer;         //!   result from shape derivatives for jet mass: 1st derivative
265   Double_t          fJetShapeMassSecondDer;        //!   result from shape derivatives for jet mass: 2nd derivative
266   Double_t          fJetShapeMassFirstSub;         //!   result from shape derivatives for jet mass: 1st order subtracted
267   Double_t          fJetShapeMassSecondSub;        //!   result from shape derivatives for jet mass: 2nd order subtracted
268   Int_t             fLabel;                        //    label to inclusive jet for constituent subtracted jet    
269
270   TArrayF           fGRNumerator;                  //!   array with angular structure function numerator
271   TArrayF           fGRDenominator;                //!   array with angular structure function denominator
272   TArrayF           fGRNumeratorSub;               //!   array with angular structure function numerator
273   TArrayF           fGRDenominatorSub;             //!   array with angular structure function denominator
274
275   Double_t          fJetShapeAngularityFirstDer;   //!   result from shape derivatives for jet Angularity: 1st derivative
276   Double_t          fJetShapeAngularitySecondDer;  //!   result from shape derivatives for jet Angularity: 2nd derivative
277   Double_t          fJetShapeAngularityFirstSub;   //!   result from shape derivatives for jet Angularity: 1st order subtracted
278   Double_t          fJetShapeAngularitySecondSub;  //!   result from shape derivatives for jet Angularity: 2nd order subtracted
279
280   Double_t          fJetShapepTDFirstDer;          //!   result from shape derivatives for jet pTD: 1st derivative
281   Double_t          fJetShapepTDSecondDer;         //!   result from shape derivatives for jet pTD: 2nd derivative
282   Double_t          fJetShapepTDFirstSub;          //!   result from shape derivatives for jet pTD: 1st order subtracted
283   Double_t          fJetShapepTDSecondSub;         //!   result from shape derivatives for jet pTD: 2nd order subtracted
284  
285   Double_t          fJetShapeCircularityFirstDer;  //!   result from shape derivatives for jet circularity: 1st derivative
286   Double_t          fJetShapeCircularitySecondDer; //!   result from shape derivatives for jet circularity: 2nd derivative
287   Double_t          fJetShapeCircularityFirstSub;  //!   result from shape derivatives for jet circularity: 1st order subtracted
288   Double_t          fJetShapeCircularitySecondSub; //!   result from shape derivatives for jetcircularity: 2nd order subtracted
289
290   Double_t          fJetShapeConstituentFirstDer;  //!   result from shape derivatives for jet const: 1st derivative
291   Double_t          fJetShapeConstituentSecondDer; //!   result from shape derivatives for jet const: 2nd derivative
292   Double_t          fJetShapeConstituentFirstSub;  //!   result from shape derivatives for jet const: 1st order subtracted
293   Double_t          fJetShapeConstituentSecondSub; //!   result from shape derivatives for jet const: 2nd order subtracted
294
295   Double_t          fJetShapeLeSubFirstDer;        //!   result from shape derivatives for jet LeSub: 1st derivative
296   Double_t          fJetShapeLeSubSecondDer;       //!   result from shape derivatives for jet LeSub: 2nd derivative
297   Double_t          fJetShapeLeSubFirstSub;        //!   result from shape derivatives for jet LeSub: 1st order subtracted
298   Double_t          fJetShapeLeSubSecondSub;       //!   result from shape derivatives for jet LeSub: 2nd order subtracted 
299
300   private:
301     struct sort_descend
302         { // sort in decreasing order
303           // first value of the pair is Pt and the second is entry index
304         bool operator () (const std::pair<Double_t, Int_t>& p1, const std::pair<Double_t, Int_t>& p2)  { return p1.first > p2.first ; }
305         };
306
307   ClassDef(AliEmcalJet,15) // Emcal jet class in cylindrical coordinates
308 };
309 #endif