ROOT logo
#ifndef AliEmcalJet_H
#define AliEmcalJet_H

#include <vector>
#include <algorithm>
#include <utility>
#include <TArrayS.h>
#include <TLorentzVector.h>
#include <TMath.h>
#include <TClonesArray.h>
#include <TVector2.h>

#include "AliVParticle.h"
#include "AliVCluster.h"
#include "AliVEvent.h"

class AliEmcalJet : public AliVParticle
{
 public:
     enum EFlavourTag{
       kDStar = 1<<0,
       kD0 = 1<<1,
       kSig1 = 1<<2,
       kSig2 = 1<<3,
       kBckgrd1 = 1<<4,
       kBckgrd2 = 1<<5,
       kBckgrd3 = 1<<6
       //.....
    };

  AliEmcalJet();
  AliEmcalJet(Double_t px, Double_t py, Double_t pz);
  AliEmcalJet(Double_t pt, Double_t eta, Double_t phi, Double_t m);
  AliEmcalJet(const AliEmcalJet &jet);
  AliEmcalJet& operator=(const AliEmcalJet &jet);

  Double_t          Px()                         const { return fPt*TMath::Cos(fPhi);  }
  Double_t          Py()                         const { return fPt*TMath::Sin(fPhi);  }
  Double_t          Pz()                         const { return fPt*TMath::SinH(fEta); }
  Double_t          Pt()                         const { return fPt;                   }
  Double_t          P()                          const { return fPt*TMath::CosH(fEta); }
  Bool_t            PxPyPz(Double_t p[3])        const { p[0]=Px();p[1]=Py();p[2]=Pz(); return 1;         }
  Double_t          Xv()                         const { return 0.;      }
  Double_t          Yv()                         const { return 0.;      }
  Double_t          Zv()                         const { return 0.;      }
  Bool_t            XvYvZv(Double_t x[3])        const { x[0]=0;x[1]=0;x[2]=0; return 1;                  }
  Double_t          OneOverPt()                  const { return 1./fPt;  }
  Double_t          Phi()                        const { return fPhi;    }
  Double_t          Theta()                      const { return 2*TMath::ATan(TMath::Exp(-fEta));         }
  Double_t          E()                          const { Double_t p=P(); return TMath::Sqrt(fM*fM+p*p); }
  Double_t          M()                          const { return fM; }
  Double_t          Eta()                        const { return fEta;    }
  Double_t          Y()                          const { Double_t e = E(); Double_t pz = Pz(); return 0.5*TMath::Log((e+pz)/(e-pz));    }
  Short_t           Charge()                     const { return 0;       }
  Int_t             GetLabel()                   const { return fLabel;  }
  Int_t             PdgCode()                    const { return 0;       }
  const Double_t   *PID()                        const { return 0;       }
  void              GetMom(TLorentzVector &vec)  const;
  void              Print(Option_t* option = "") const;

  Double_t          Area()                       const { return fArea;                     }
  Double_t          AreaPt()                     const { return fArea;                     }
  Double_t          AreaEta()                    const { return fAreaEta;                  }
  Double_t          AreaPhi()                    const { return fAreaPhi;                  }
  Double_t          AreaE()                      const { return fAreaE;                    }
  Double_t          AreaEmc()                    const { return fAreaEmc;                  }
  Bool_t            AxisInEmcal()                const { return fAxisInEmcal;              }
  Int_t             Compare(const TObject* obj)  const;
  Short_t           ClusterAt(Int_t idx)         const { return fClusterIDs.At(idx);       }
  AliVCluster      *ClusterAt(Int_t idx, TClonesArray *ca)  const { if (!ca) return 0; return dynamic_cast<AliVCluster*>(ca->At(ClusterAt(idx))); }
  AliVCluster      *GetLeadingCluster(TClonesArray *clusters) const;
  UShort_t          GetNumberOfClusters()        const { return fClusterIDs.GetSize();     }
  UShort_t          GetNumberOfTracks()          const { return fTrackIDs.GetSize();       }
  UShort_t          GetNumberOfConstituents()    const { return GetNumberOfClusters()+GetNumberOfTracks();       }
  Double_t          FracEmcalArea()              const { return fAreaEmc/fArea;            }
  Bool_t            IsInsideEmcal()              const { return (fAreaEmc/fArea>0.999);    }
  Bool_t            IsInEmcal()                  const { return (fAreaEmc>0);              }
  Bool_t            IsMC()                       const { return (Bool_t)(MCPt() > 0);      }
  Bool_t            IsSortable()                 const { return kTRUE;                     }
  Double_t          MaxNeutralPt()               const { return fMaxNPt;                   }
  Double_t          MaxChargedPt()               const { return fMaxCPt;                   }
  Double_t          NEF()                        const { return fNEF;                      }
  UShort_t          Nn()                         const { return fNn;                       }
  UShort_t          Nch()                        const { return fNch;                      }
  UShort_t          N()                          const { return Nch()+Nn();                }
  Int_t             NEmc()                       const { return fNEmc;                     }
  Double_t          MCPt()                       const { return fMCPt;                     }
  Double_t          MaxClusterPt()               const { return MaxNeutralPt();            }
  Double_t          MaxTrackPt()                 const { return MaxChargedPt();            }
  Double_t          MaxPartPt()                  const { return fMaxCPt < fMaxNPt ? fMaxNPt : fMaxCPt;     }
  Double_t          PtEmc()                      const { return fPtEmc;                    }
  Double_t          PtSub()                      const { return fPtSub;                    }
  Double_t          PtSubVect()                  const { return fPtSubVect;                }
  Double_t          PtSub(Double_t rho, Bool_t save = kFALSE);
  Double_t          PtSubVect(Double_t rho, Bool_t save = kFALSE);
  TLorentzVector    SubtractRhoVect(Double_t rho, Bool_t save = kFALSE);
  Short_t           TrackAt(Int_t idx)           const { return fTrackIDs.At(idx);         }
  AliVParticle     *TrackAt(Int_t idx, TClonesArray *ta)  const { if (!ta) return 0; return dynamic_cast<AliVParticle*>(ta->At(TrackAt(idx))); }
  AliVParticle     *GetLeadingTrack(TClonesArray *tracks) const;
  Int_t             GetFlavour()                 const { return fFlavourTagging;           }

  void              AddClusterAt(Int_t clus, Int_t idx){ fClusterIDs.AddAt(clus, idx);     }
  void              AddFlavourTag(Int_t tag)           { fFlavourTagging |= tag; }
  void              AddTrackAt(Int_t track, Int_t idx) { fTrackIDs.AddAt(track, idx);      }
  void              Clear(Option_t */*option*/="")     { fClusterIDs.Set(0); fTrackIDs.Set(0); fClosestJets[0] = 0; fClosestJets[1] = 0;
                                                         fClosestJetsDist[0] = 0; fClosestJetsDist[1] = 0; fMatched = 0; fPtSub = 0; }
  Double_t          DeltaR(const AliVParticle* part) const;
  Double_t          GetZ ( const Double_t trkPx, const Double_t trkPy, const Double_t trkPz ) const; // Get Z of constituent trk
  Double_t          GetZ ( const AliVParticle* trk )       const; // Get Z of constituent trk
  Double_t          GetXi ( const AliVParticle* trk )      const { return TMath::Log ( 1/GetZ (trk) ); } // Get Xi of constituent trk
  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

  void              SetLabel(Int_t l)                  { fLabel = l;                       }
  void              SetArea(Double_t a)                { fArea    = a;                     }
  void              SetAreaEta(Double_t a)             { fAreaEta = a;                     }
  void              SetAreaPhi(Double_t a)             { fAreaPhi = TVector2::Phi_0_2pi(a); }
  void              SetAreaE(Double_t a)               { fAreaE = a;                       }
  void              SetAreaEmc(Double_t a)             { fAreaEmc = a;                     }
  void              SetAxisInEmcal(Bool_t b)           { fAxisInEmcal = b;                 }
  void              SetFlavour(Int_t flavour)          { fFlavourTagging = flavour;        }
  void              SetMaxNeutralPt(Double32_t t)      { fMaxNPt  = t;                     }
  void              SetMaxChargedPt(Double32_t t)      { fMaxCPt  = t;                     }
  void              SetNEF(Double_t nef)               { fNEF     = nef;                   }
  void              SetNumberOfClusters(Int_t n)       { fClusterIDs.Set(n);               }
  void              SetNumberOfTracks(Int_t n)         { fTrackIDs.Set(n);                 }
  void              SetNumberOfCharged(Int_t n)        { fNch = n;                         }
  void              SetNumberOfNeutrals(Int_t n)       { fNn = n;                          }
  void              SetMCPt(Double_t p)                { fMCPt = p;                        }
  void              SortConstituents();
  std::vector<int>  SortConstituentsPt(TClonesArray *tracks) const;
  void              SetNEmc(Int_t n)                   { fNEmc           = n;              }
  void              SetPtEmc(Double_t pt)              { fPtEmc          = pt;             }
  void              SetPtSub(Double_t ps)              { fPtSub          = ps;             }
  void              SetPtSubVect(Double_t ps)          { fPtSubVect      = ps;             }
  Bool_t            TestFlavourTag(Int_t tag)    const { return (Bool_t)((tag & fFlavourTagging) !=0); }

  // Trigger
  Bool_t            IsTriggerJet(UInt_t trigger=AliVEvent::kEMCEJE) const   { return (Bool_t)((fTriggers & trigger) != 0); }
  void              SetTrigger(UInt_t trigger)                              { fTriggers  = trigger;                        }
  void              AddTrigger(UInt_t trigger)                              { fTriggers |= trigger;                        }

  // Matching
  void              SetClosestJet(AliEmcalJet *j, Double_t d)       { fClosestJets[0] = j; fClosestJetsDist[0] = d    ; }
  void              SetSecondClosestJet(AliEmcalJet *j, Double_t d) { fClosestJets[1] = j; fClosestJetsDist[1] = d    ; }
  void              SetMatchedToClosest(UShort_t m)                 { fMatched        = 0; fMatchingType       = m    ; }
  void              SetMatchedToSecondClosest(UShort_t m)           { fMatched        = 1; fMatchingType       = m    ; }
  void              ResetMatching();
  AliEmcalJet*      ClosestJet()                              const { return fClosestJets[0]                          ; }
  Double_t          ClosestJetDistance()                      const { return fClosestJetsDist[0]                      ; }
  AliEmcalJet*      SecondClosestJet()                        const { return fClosestJets[1]                          ; }
  Double_t          SecondClosestJetDistance()                const { return fClosestJetsDist[1]                      ; }
  AliEmcalJet*      MatchedJet()                              const { return fMatched < 2 ? fClosestJets[fMatched] : 0; }
  UShort_t          GetMatchingType()                         const { return fMatchingType                            ; }

  void              SetTaggedJet(AliEmcalJet *j)                    { fTaggedJet = j                                  ; }
  void              SetTagStatus(Int_t i)                           { fTagStatus = i                                  ; }
  AliEmcalJet*      GetTaggedJet()                            const { return fTaggedJet                               ; }
  Int_t             GetTagStatus()                            const { return fTagStatus                               ; }

  //jet shape derivatives
  //jet mass
  void              SetFirstDerivative(Double_t d)                  { fJetShapeMassFirstDer = d                       ; }
  void              SetSecondDerivative(Double_t d)                 { fJetShapeMassSecondDer = d                      ; }
  void              SetFirstOrderSubtracted(Double_t d)             { fJetShapeMassFirstSub = d                       ; }
  void              SetSecondOrderSubtracted(Double_t d)            { fJetShapeMassSecondSub = d                      ; }
  Double_t          GetFirstDerivative()                      const { return fJetShapeMassFirstDer                    ; }
  Double_t          GetSecondDerivative()                     const { return fJetShapeMassSecondDer                   ; }
  Double_t          GetFirstOrderSubtracted()                 const { return fJetShapeMassFirstSub                    ; }
  Double_t          GetSecondOrderSubtracted()                const { return fJetShapeMassSecondSub                   ; }

  //jet structure function
  TArrayF           GetGRNumerator()                          const { return fGRNumerator                             ; }
  TArrayF           GetGRDenominator()                        const { return fGRDenominator                           ; }
  TArrayF           GetGRNumeratorSub()                       const { return fGRNumeratorSub                          ; }
  TArrayF           GetGRDenominatorSub()                     const { return fGRDenominatorSub                        ; }
  void              AddGRNumAt(Float_t num, Int_t idx)              { fGRNumerator.AddAt(num, idx)                    ; }
  void              AddGRDenAt(Float_t den, Int_t idx)              { fGRDenominator.AddAt(den, idx)                  ; }
  void              SetGRNumSize(UInt_t s)                          { fGRNumerator.Set(s)                             ; }
  void              SetGRDenSize(UInt_t s)                          { fGRDenominator.Set(s)                           ; }

  void              AddGRNumSubAt(Float_t num, Int_t idx)           { fGRNumeratorSub.AddAt(num, idx)                 ; }
  void              AddGRDenSubAt(Float_t den, Int_t idx)           { fGRDenominatorSub.AddAt(den, idx)               ; }
  void              SetGRNumSubSize(UInt_t s)                       { fGRNumeratorSub.Set(s)                          ; }
  void              SetGRDenSubSize(UInt_t s)                       { fGRDenominatorSub.Set(s)                        ; }
  void              PrintGR();

  //Angularity
  void              SetFirstDerivativeAngularity(Double_t d)        { fJetShapeAngularityFirstDer = d                 ; }
  void              SetSecondDerivativeAngularity(Double_t d)       { fJetShapeAngularitySecondDer = d                ; }
  void              SetFirstOrderSubtractedAngularity(Double_t d)   { fJetShapeAngularityFirstSub = d                 ; }
  void              SetSecondOrderSubtractedAngularity(Double_t d)  { fJetShapeAngularitySecondSub = d                ; }
  Double_t          GetFirstDerivativeAngularity()            const { return fJetShapeAngularityFirstDer              ; }
  Double_t          GetSecondDerivativeAngularity()           const { return fJetShapeAngularitySecondDer             ; }
  Double_t          GetFirstOrderSubtractedAngularity()       const { return fJetShapeAngularityFirstSub              ; }
  Double_t          GetSecondOrderSubtractedAngularity()      const { return fJetShapeAngularitySecondSub             ; }

  //pTD
  void              SetFirstDerivativepTD(Double_t d)               { fJetShapepTDFirstDer = d                        ; }
  void              SetSecondDerivativepTD(Double_t d)              { fJetShapepTDSecondDer = d                       ; }
  void              SetFirstOrderSubtractedpTD(Double_t d)          { fJetShapepTDFirstSub = d                        ; }
  void              SetSecondOrderSubtractedpTD(Double_t d)         { fJetShapepTDSecondSub = d                       ; }
  Double_t          GetFirstDerivativepTD()                   const { return fJetShapepTDFirstDer                     ; }
  Double_t          GetSecondDerivativepTD()                  const { return fJetShapepTDSecondDer                    ; }
  Double_t          GetFirstOrderSubtractedpTD()              const { return fJetShapepTDFirstSub                     ; }
  Double_t          GetSecondOrderSubtractedpTD()             const { return fJetShapepTDSecondSub                    ; }

  //Circularity
  void              SetFirstDerivativeCircularity(Double_t d)       { fJetShapeCircularityFirstDer = d                ; }
  void              SetSecondDerivativeCircularity(Double_t d)      { fJetShapeCircularitySecondDer = d               ; }
  void              SetFirstOrderSubtractedCircularity(Double_t d)  { fJetShapeCircularityFirstSub = d                ; }
  void              SetSecondOrderSubtractedCircularity(Double_t d) { fJetShapeCircularitySecondSub = d               ; }
  Double_t          GetFirstDerivativeCircularity()           const { return fJetShapeCircularityFirstDer             ; }
  Double_t          GetSecondDerivativeCircularity()          const { return fJetShapeCircularitySecondDer            ; }
  Double_t          GetFirstOrderSubtractedCircularity()      const { return fJetShapeCircularityFirstSub             ; }
  Double_t          GetSecondOrderSubtractedCircularity()     const { return fJetShapeCircularitySecondSub            ; }

  //Sigma2
  void              SetFirstDerivativeSigma2(Double_t d)            { fJetShapeSigma2FirstDer = d                     ; }
  void              SetSecondDerivativeSigma2(Double_t d)           { fJetShapeSigma2SecondDer = d                    ; }
  void              SetFirstOrderSubtractedSigma2(Double_t d)       { fJetShapeSigma2FirstSub = d                     ; }
  void              SetSecondOrderSubtractedSigma2(Double_t d)      { fJetShapeSigma2SecondSub = d                    ; }
  Double_t          GetFirstDerivativeSigma2()           const      { return fJetShapeSigma2FirstDer                  ; }
  Double_t          GetSecondDerivativeSigma2()          const      { return fJetShapeSigma2SecondDer                 ; }
  Double_t          GetFirstOrderSubtractedSigma2()      const      { return fJetShapeSigma2FirstSub                  ; }
  Double_t          GetSecondOrderSubtractedSigma2()     const      { return fJetShapeSigma2SecondSub                 ; }


  //number of contituents
  void              SetFirstDerivativeConstituent(Double_t d)       { fJetShapeConstituentFirstDer = d                ; }
  void              SetSecondDerivativeConstituent(Double_t d)      { fJetShapeConstituentSecondDer = d               ; }
  void              SetFirstOrderSubtractedConstituent(Double_t d)  { fJetShapeConstituentFirstSub = d                ; }
  void              SetSecondOrderSubtractedConstituent(Double_t d) { fJetShapeConstituentSecondSub = d               ; }
  Double_t          GetFirstDerivativeConstituent()           const { return fJetShapeConstituentFirstDer             ; }
  Double_t          GetSecondDerivativeConstituent()          const { return fJetShapeConstituentSecondDer            ; }
  Double_t          GetFirstOrderSubtractedConstituent()      const { return fJetShapeConstituentFirstSub             ; }
  Double_t          GetSecondOrderSubtractedConstituent()     const { return fJetShapeConstituentSecondSub            ; }

  //leading minus subleading constituent
  void              SetFirstDerivativeLeSub(Double_t d)             { fJetShapeLeSubFirstDer = d                      ; }
  void              SetSecondDerivativeLeSub(Double_t d)            { fJetShapeLeSubSecondDer = d                     ; }
  void              SetFirstOrderSubtractedLeSub(Double_t d)        { fJetShapeLeSubFirstSub = d                      ; }
  void              SetSecondOrderSubtractedLeSub(Double_t d)       { fJetShapeLeSubSecondSub = d                     ; }
  Double_t          GetFirstDerivativeLeSub()                 const { return fJetShapeLeSubFirstDer                   ; }
  Double_t          GetSecondDerivativeLeSub()                const { return fJetShapeLeSubSecondDer                  ; }
  Double_t          GetFirstOrderSubtractedLeSub()            const { return fJetShapeLeSubFirstSub                   ; }
  Double_t          GetSecondOrderSubtractedLeSub()           const { return fJetShapeLeSubSecondSub                  ; }

 protected:
  Double32_t        fPt;                  //[0,0,12]   pt
  Double32_t        fEta;                 //[-1,1,12]  eta
  Double32_t        fPhi;                 //[0,6.3,12] phi
  Double32_t        fM;                   //[0,0,8]    mass
  Double32_t        fNEF;                 //[0,1,8]    neutral energy fraction
  Double32_t        fArea;                //[0,0,12]   area (transverse)
  Double32_t        fAreaEta;             //[0,0,12]   area eta
  Double32_t        fAreaPhi;             //[0,0,12]   area phi
  Double32_t        fAreaE;               //[0,0,12]   temporal area component
  Double32_t        fAreaEmc;             //[0,0,12]   area on EMCAL surface (determined from ghosts)
  Bool_t            fAxisInEmcal;         //           =true if jet axis inside EMCAL acceptance
  Int_t             fFlavourTagging;      // tag jet with a falvour, bit 0 = no tag; bit 1= Dstar; bit 2 = D0
  Double32_t        fMaxCPt;              //[0,0,12]   pt of maximum charged constituent
  Double32_t        fMaxNPt;              //[0,0,12]   pt of maximum neutral constituent
  Double32_t        fMCPt;                //           pt from MC particles contributing to the jet
  Int_t             fNn;                  //           number of neutral constituents
  Int_t             fNch;                 //           number of charged constituents
  Double32_t        fPtEmc;               //[0,0,12]   pt in EMCAL acceptance
  Int_t             fNEmc;                //           number of constituents in EMCAL acceptance
  TArrayI           fClusterIDs;          //           array containing ids of cluster constituents
  TArrayI           fTrackIDs;            //           array containing ids of track constituents
  AliEmcalJet      *fClosestJets[2];      //!          if this is MC it contains the two closest detector level jets in order of distance and viceversa
  Double32_t        fClosestJetsDist[2];  //!          distance to closest jets (see above)
  UShort_t          fMatched;             //!          0,1 if it is matched with one of the closest jets; 2 if it is not matched
  UShort_t          fMatchingType;        //!          matching type
  AliEmcalJet      *fTaggedJet;           //!          jet tagged to this jet
  Int_t             fTagStatus;           //!          status of tagging -1: NA 0: not tagged 1: tagged
  Double_t          fPtSub;               //!          background subtracted pt (not stored set from outside)
  Double_t          fPtSubVect;           //!          background vector subtracted pt (not stored set from outside)
  UInt_t            fTriggers;            //!          triggers that the jet might have fired (AliVEvent::EOfflineTriggerTypes)

  Double_t          fJetShapeMassFirstDer;         //!   result from shape derivatives for jet mass: 1st derivative
  Double_t          fJetShapeMassSecondDer;        //!   result from shape derivatives for jet mass: 2nd derivative
  Double_t          fJetShapeMassFirstSub;         //!   result from shape derivatives for jet mass: 1st order subtracted
  Double_t          fJetShapeMassSecondSub;        //!   result from shape derivatives for jet mass: 2nd order subtracted
  Int_t             fLabel;                        //    label to inclusive jet for constituent subtracted jet

  TArrayF           fGRNumerator;                  //!   array with angular structure function numerator
  TArrayF           fGRDenominator;                //!   array with angular structure function denominator
  TArrayF           fGRNumeratorSub;               //!   array with angular structure function numerator
  TArrayF           fGRDenominatorSub;             //!   array with angular structure function denominator

  Double_t          fJetShapeAngularityFirstDer;   //!   result from shape derivatives for jet Angularity: 1st derivative
  Double_t          fJetShapeAngularitySecondDer;  //!   result from shape derivatives for jet Angularity: 2nd derivative
  Double_t          fJetShapeAngularityFirstSub;   //!   result from shape derivatives for jet Angularity: 1st order subtracted
  Double_t          fJetShapeAngularitySecondSub;  //!   result from shape derivatives for jet Angularity: 2nd order subtracted

  Double_t          fJetShapepTDFirstDer;          //!   result from shape derivatives for jet pTD: 1st derivative
  Double_t          fJetShapepTDSecondDer;         //!   result from shape derivatives for jet pTD: 2nd derivative
  Double_t          fJetShapepTDFirstSub;          //!   result from shape derivatives for jet pTD: 1st order subtracted
  Double_t          fJetShapepTDSecondSub;         //!   result from shape derivatives for jet pTD: 2nd order subtracted

  Double_t          fJetShapeCircularityFirstDer;  //!   result from shape derivatives for jet circularity: 1st derivative
  Double_t          fJetShapeCircularitySecondDer; //!   result from shape derivatives for jet circularity: 2nd derivative
  Double_t          fJetShapeCircularityFirstSub;  //!   result from shape derivatives for jet circularity: 1st order subtracted
  Double_t          fJetShapeCircularitySecondSub; //!   result from shape derivatives for jetcircularity: 2nd order subtracted

  Double_t          fJetShapeSigma2FirstDer;       //!   result from shape derivatives for jet sigma2: 1st derivative
  Double_t          fJetShapeSigma2SecondDer;      //!   result from shape derivatives for jet sigma2: 2nd derivative
  Double_t          fJetShapeSigma2FirstSub;       //!   result from shape derivatives for jet sigma2: 1st order subtracted
  Double_t          fJetShapeSigma2SecondSub;      //!   result from shape derivatives for jetsigma2: 2nd order subtracted

  Double_t          fJetShapeConstituentFirstDer;  //!   result from shape derivatives for jet const: 1st derivative
  Double_t          fJetShapeConstituentSecondDer; //!   result from shape derivatives for jet const: 2nd derivative
  Double_t          fJetShapeConstituentFirstSub;  //!   result from shape derivatives for jet const: 1st order subtracted
  Double_t          fJetShapeConstituentSecondSub; //!   result from shape derivatives for jet const: 2nd order subtracted

  Double_t          fJetShapeLeSubFirstDer;        //!   result from shape derivatives for jet LeSub: 1st derivative
  Double_t          fJetShapeLeSubSecondDer;       //!   result from shape derivatives for jet LeSub: 2nd derivative
  Double_t          fJetShapeLeSubFirstSub;        //!   result from shape derivatives for jet LeSub: 1st order subtracted
  Double_t          fJetShapeLeSubSecondSub;       //!   result from shape derivatives for jet LeSub: 2nd order subtracted

  private:
    struct sort_descend
        { // sort in decreasing order
          // first value of the pair is Pt and the second is entry index
        bool operator () (const std::pair<Double_t, Int_t>& p1, const std::pair<Double_t, Int_t>& p2)  { return p1.first > p2.first ; }
        };

  ClassDef(AliEmcalJet,16) // Emcal jet class in cylindrical coordinates
};
#endif
 AliEmcalJet.h:1
 AliEmcalJet.h:2
 AliEmcalJet.h:3
 AliEmcalJet.h:4
 AliEmcalJet.h:5
 AliEmcalJet.h:6
 AliEmcalJet.h:7
 AliEmcalJet.h:8
 AliEmcalJet.h:9
 AliEmcalJet.h:10
 AliEmcalJet.h:11
 AliEmcalJet.h:12
 AliEmcalJet.h:13
 AliEmcalJet.h:14
 AliEmcalJet.h:15
 AliEmcalJet.h:16
 AliEmcalJet.h:17
 AliEmcalJet.h:18
 AliEmcalJet.h:19
 AliEmcalJet.h:20
 AliEmcalJet.h:21
 AliEmcalJet.h:22
 AliEmcalJet.h:23
 AliEmcalJet.h:24
 AliEmcalJet.h:25
 AliEmcalJet.h:26
 AliEmcalJet.h:27
 AliEmcalJet.h:28
 AliEmcalJet.h:29
 AliEmcalJet.h:30
 AliEmcalJet.h:31
 AliEmcalJet.h:32
 AliEmcalJet.h:33
 AliEmcalJet.h:34
 AliEmcalJet.h:35
 AliEmcalJet.h:36
 AliEmcalJet.h:37
 AliEmcalJet.h:38
 AliEmcalJet.h:39
 AliEmcalJet.h:40
 AliEmcalJet.h:41
 AliEmcalJet.h:42
 AliEmcalJet.h:43
 AliEmcalJet.h:44
 AliEmcalJet.h:45
 AliEmcalJet.h:46
 AliEmcalJet.h:47
 AliEmcalJet.h:48
 AliEmcalJet.h:49
 AliEmcalJet.h:50
 AliEmcalJet.h:51
 AliEmcalJet.h:52
 AliEmcalJet.h:53
 AliEmcalJet.h:54
 AliEmcalJet.h:55
 AliEmcalJet.h:56
 AliEmcalJet.h:57
 AliEmcalJet.h:58
 AliEmcalJet.h:59
 AliEmcalJet.h:60
 AliEmcalJet.h:61
 AliEmcalJet.h:62
 AliEmcalJet.h:63
 AliEmcalJet.h:64
 AliEmcalJet.h:65
 AliEmcalJet.h:66
 AliEmcalJet.h:67
 AliEmcalJet.h:68
 AliEmcalJet.h:69
 AliEmcalJet.h:70
 AliEmcalJet.h:71
 AliEmcalJet.h:72
 AliEmcalJet.h:73
 AliEmcalJet.h:74
 AliEmcalJet.h:75
 AliEmcalJet.h:76
 AliEmcalJet.h:77
 AliEmcalJet.h:78
 AliEmcalJet.h:79
 AliEmcalJet.h:80
 AliEmcalJet.h:81
 AliEmcalJet.h:82
 AliEmcalJet.h:83
 AliEmcalJet.h:84
 AliEmcalJet.h:85
 AliEmcalJet.h:86
 AliEmcalJet.h:87
 AliEmcalJet.h:88
 AliEmcalJet.h:89
 AliEmcalJet.h:90
 AliEmcalJet.h:91
 AliEmcalJet.h:92
 AliEmcalJet.h:93
 AliEmcalJet.h:94
 AliEmcalJet.h:95
 AliEmcalJet.h:96
 AliEmcalJet.h:97
 AliEmcalJet.h:98
 AliEmcalJet.h:99
 AliEmcalJet.h:100
 AliEmcalJet.h:101
 AliEmcalJet.h:102
 AliEmcalJet.h:103
 AliEmcalJet.h:104
 AliEmcalJet.h:105
 AliEmcalJet.h:106
 AliEmcalJet.h:107
 AliEmcalJet.h:108
 AliEmcalJet.h:109
 AliEmcalJet.h:110
 AliEmcalJet.h:111
 AliEmcalJet.h:112
 AliEmcalJet.h:113
 AliEmcalJet.h:114
 AliEmcalJet.h:115
 AliEmcalJet.h:116
 AliEmcalJet.h:117
 AliEmcalJet.h:118
 AliEmcalJet.h:119
 AliEmcalJet.h:120
 AliEmcalJet.h:121
 AliEmcalJet.h:122
 AliEmcalJet.h:123
 AliEmcalJet.h:124
 AliEmcalJet.h:125
 AliEmcalJet.h:126
 AliEmcalJet.h:127
 AliEmcalJet.h:128
 AliEmcalJet.h:129
 AliEmcalJet.h:130
 AliEmcalJet.h:131
 AliEmcalJet.h:132
 AliEmcalJet.h:133
 AliEmcalJet.h:134
 AliEmcalJet.h:135
 AliEmcalJet.h:136
 AliEmcalJet.h:137
 AliEmcalJet.h:138
 AliEmcalJet.h:139
 AliEmcalJet.h:140
 AliEmcalJet.h:141
 AliEmcalJet.h:142
 AliEmcalJet.h:143
 AliEmcalJet.h:144
 AliEmcalJet.h:145
 AliEmcalJet.h:146
 AliEmcalJet.h:147
 AliEmcalJet.h:148
 AliEmcalJet.h:149
 AliEmcalJet.h:150
 AliEmcalJet.h:151
 AliEmcalJet.h:152
 AliEmcalJet.h:153
 AliEmcalJet.h:154
 AliEmcalJet.h:155
 AliEmcalJet.h:156
 AliEmcalJet.h:157
 AliEmcalJet.h:158
 AliEmcalJet.h:159
 AliEmcalJet.h:160
 AliEmcalJet.h:161
 AliEmcalJet.h:162
 AliEmcalJet.h:163
 AliEmcalJet.h:164
 AliEmcalJet.h:165
 AliEmcalJet.h:166
 AliEmcalJet.h:167
 AliEmcalJet.h:168
 AliEmcalJet.h:169
 AliEmcalJet.h:170
 AliEmcalJet.h:171
 AliEmcalJet.h:172
 AliEmcalJet.h:173
 AliEmcalJet.h:174
 AliEmcalJet.h:175
 AliEmcalJet.h:176
 AliEmcalJet.h:177
 AliEmcalJet.h:178
 AliEmcalJet.h:179
 AliEmcalJet.h:180
 AliEmcalJet.h:181
 AliEmcalJet.h:182
 AliEmcalJet.h:183
 AliEmcalJet.h:184
 AliEmcalJet.h:185
 AliEmcalJet.h:186
 AliEmcalJet.h:187
 AliEmcalJet.h:188
 AliEmcalJet.h:189
 AliEmcalJet.h:190
 AliEmcalJet.h:191
 AliEmcalJet.h:192
 AliEmcalJet.h:193
 AliEmcalJet.h:194
 AliEmcalJet.h:195
 AliEmcalJet.h:196
 AliEmcalJet.h:197
 AliEmcalJet.h:198
 AliEmcalJet.h:199
 AliEmcalJet.h:200
 AliEmcalJet.h:201
 AliEmcalJet.h:202
 AliEmcalJet.h:203
 AliEmcalJet.h:204
 AliEmcalJet.h:205
 AliEmcalJet.h:206
 AliEmcalJet.h:207
 AliEmcalJet.h:208
 AliEmcalJet.h:209
 AliEmcalJet.h:210
 AliEmcalJet.h:211
 AliEmcalJet.h:212
 AliEmcalJet.h:213
 AliEmcalJet.h:214
 AliEmcalJet.h:215
 AliEmcalJet.h:216
 AliEmcalJet.h:217
 AliEmcalJet.h:218
 AliEmcalJet.h:219
 AliEmcalJet.h:220
 AliEmcalJet.h:221
 AliEmcalJet.h:222
 AliEmcalJet.h:223
 AliEmcalJet.h:224
 AliEmcalJet.h:225
 AliEmcalJet.h:226
 AliEmcalJet.h:227
 AliEmcalJet.h:228
 AliEmcalJet.h:229
 AliEmcalJet.h:230
 AliEmcalJet.h:231
 AliEmcalJet.h:232
 AliEmcalJet.h:233
 AliEmcalJet.h:234
 AliEmcalJet.h:235
 AliEmcalJet.h:236
 AliEmcalJet.h:237
 AliEmcalJet.h:238
 AliEmcalJet.h:239
 AliEmcalJet.h:240
 AliEmcalJet.h:241
 AliEmcalJet.h:242
 AliEmcalJet.h:243
 AliEmcalJet.h:244
 AliEmcalJet.h:245
 AliEmcalJet.h:246
 AliEmcalJet.h:247
 AliEmcalJet.h:248
 AliEmcalJet.h:249
 AliEmcalJet.h:250
 AliEmcalJet.h:251
 AliEmcalJet.h:252
 AliEmcalJet.h:253
 AliEmcalJet.h:254
 AliEmcalJet.h:255
 AliEmcalJet.h:256
 AliEmcalJet.h:257
 AliEmcalJet.h:258
 AliEmcalJet.h:259
 AliEmcalJet.h:260
 AliEmcalJet.h:261
 AliEmcalJet.h:262
 AliEmcalJet.h:263
 AliEmcalJet.h:264
 AliEmcalJet.h:265
 AliEmcalJet.h:266
 AliEmcalJet.h:267
 AliEmcalJet.h:268
 AliEmcalJet.h:269
 AliEmcalJet.h:270
 AliEmcalJet.h:271
 AliEmcalJet.h:272
 AliEmcalJet.h:273
 AliEmcalJet.h:274
 AliEmcalJet.h:275
 AliEmcalJet.h:276
 AliEmcalJet.h:277
 AliEmcalJet.h:278
 AliEmcalJet.h:279
 AliEmcalJet.h:280
 AliEmcalJet.h:281
 AliEmcalJet.h:282
 AliEmcalJet.h:283
 AliEmcalJet.h:284
 AliEmcalJet.h:285
 AliEmcalJet.h:286
 AliEmcalJet.h:287
 AliEmcalJet.h:288
 AliEmcalJet.h:289
 AliEmcalJet.h:290
 AliEmcalJet.h:291
 AliEmcalJet.h:292
 AliEmcalJet.h:293
 AliEmcalJet.h:294
 AliEmcalJet.h:295
 AliEmcalJet.h:296
 AliEmcalJet.h:297
 AliEmcalJet.h:298
 AliEmcalJet.h:299
 AliEmcalJet.h:300
 AliEmcalJet.h:301
 AliEmcalJet.h:302
 AliEmcalJet.h:303
 AliEmcalJet.h:304
 AliEmcalJet.h:305
 AliEmcalJet.h:306
 AliEmcalJet.h:307
 AliEmcalJet.h:308
 AliEmcalJet.h:309
 AliEmcalJet.h:310
 AliEmcalJet.h:311
 AliEmcalJet.h:312
 AliEmcalJet.h:313
 AliEmcalJet.h:314
 AliEmcalJet.h:315
 AliEmcalJet.h:316
 AliEmcalJet.h:317
 AliEmcalJet.h:318
 AliEmcalJet.h:319
 AliEmcalJet.h:320
 AliEmcalJet.h:321
 AliEmcalJet.h:322
 AliEmcalJet.h:323
 AliEmcalJet.h:324
 AliEmcalJet.h:325
 AliEmcalJet.h:326
 AliEmcalJet.h:327
 AliEmcalJet.h:328
 AliEmcalJet.h:329
 AliEmcalJet.h:330