#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 *="") { 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;
Double_t GetZ ( const AliVParticle* trk ) const;
Double_t GetXi ( const AliVParticle* trk ) const { return TMath::Log ( 1/GetZ (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 ) ); }
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); }
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; }
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 ; }
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 ; }
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();
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 ; }
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 ; }
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 ; }
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 ; }
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 ; }
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;
Double32_t fEta;
Double32_t fPhi;
Double32_t fM;
Double32_t fNEF;
Double32_t fArea;
Double32_t fAreaEta;
Double32_t fAreaPhi;
Double32_t fAreaE;
Double32_t fAreaEmc;
Bool_t fAxisInEmcal;
Int_t fFlavourTagging;
Double32_t fMaxCPt;
Double32_t fMaxNPt;
Double32_t fMCPt;
Int_t fNn;
Int_t fNch;
Double32_t fPtEmc;
Int_t fNEmc;
TArrayI fClusterIDs;
TArrayI fTrackIDs;
AliEmcalJet *fClosestJets[2];
Double32_t fClosestJetsDist[2];
UShort_t fMatched;
UShort_t fMatchingType;
AliEmcalJet *fTaggedJet;
Int_t fTagStatus;
Double_t fPtSub;
Double_t fPtSubVect;
UInt_t fTriggers;
Double_t fJetShapeMassFirstDer;
Double_t fJetShapeMassSecondDer;
Double_t fJetShapeMassFirstSub;
Double_t fJetShapeMassSecondSub;
Int_t fLabel;
TArrayF fGRNumerator;
TArrayF fGRDenominator;
TArrayF fGRNumeratorSub;
TArrayF fGRDenominatorSub;
Double_t fJetShapeAngularityFirstDer;
Double_t fJetShapeAngularitySecondDer;
Double_t fJetShapeAngularityFirstSub;
Double_t fJetShapeAngularitySecondSub;
Double_t fJetShapepTDFirstDer;
Double_t fJetShapepTDSecondDer;
Double_t fJetShapepTDFirstSub;
Double_t fJetShapepTDSecondSub;
Double_t fJetShapeCircularityFirstDer;
Double_t fJetShapeCircularitySecondDer;
Double_t fJetShapeCircularityFirstSub;
Double_t fJetShapeCircularitySecondSub;
Double_t fJetShapeSigma2FirstDer;
Double_t fJetShapeSigma2SecondDer;
Double_t fJetShapeSigma2FirstSub;
Double_t fJetShapeSigma2SecondSub;
Double_t fJetShapeConstituentFirstDer;
Double_t fJetShapeConstituentSecondDer;
Double_t fJetShapeConstituentFirstSub;
Double_t fJetShapeConstituentSecondSub;
Double_t fJetShapeLeSubFirstDer;
Double_t fJetShapeLeSubSecondDer;
Double_t fJetShapeLeSubFirstSub;
Double_t fJetShapeLeSubSecondSub;
private:
struct sort_descend
{
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)
};
#endif