//  
// Author: Oliver Maria Kind <mailto: kind@mail.desy.de>
// Update: $Id: HepJet.h,v 1.16 2016/12/22 20:38:20 kind Exp $
// Copyright: 2008 (C) Oliver Maria Kind
//
#ifndef HEP_HepJet
#define HEP_HepJet
#ifndef ROOT_TObject
#include <TObject.h>
#endif
#ifndef ROOT_TLorentzVector
#include <TLorentzVector.h>
#endif

class HepParticle;

class HepJet : public TObject {
    
  protected:
    Int_t           fId;  // Id number (for convenience only)
    TLorentzVector  fP;   // Jet 4-momentum (GeV)
    
  private:
    Bool_t          fCompute;  //! Compute the transient variables at the next possibility
    Float_t         fEt;       //! Transverse energy (GeV)
    Float_t         fPt;       //! Transverse momentum (GeV)
    Float_t         fTheta;    //! Polar angle of the 3-momentum (rad)
    Float_t         fCosTheta; //! Cosine of the polar angle
    Float_t         fPhi;      //! Azimuth of the 3-momentum (rad)
    Float_t         fEta;      //! Pseudo-rapidity of the 3-momentum
 
  public:
    HepJet();
    HepJet(Int_t Id, Float_t E, Float_t Px, Float_t Py, Float_t Pz);
    virtual ~HepJet();
    virtual void Clear(Option_t *option = "");
    virtual Int_t Compare(const TObject *obj) const;
    virtual void Print(Option_t *option = "");
    static void PrintHeader();
    static void PrintFooter();
    Float_t DeltaPhi(HepParticle *prt) const;
    Float_t DeltaPhi(HepJet *jet) const;
    Float_t DeltaEta(HepParticle *prt);
    Float_t DeltaEta(HepJet *jet);
    Float_t DeltaR(HepParticle *prt);
    Float_t DeltaR(HepJet *jet);
    Float_t DeltaPtFrac(HepParticle *prt);
    Float_t DeltaPtFrac(HepJet *jet);
    Float_t DeltaEtFrac(HepParticle *prt);
    Float_t DeltaEtFrac(HepJet *jet);
    Float_t InvDeltaPtFrac(HepParticle *prt);
    Double_t PtRel(HepParticle* prt);
    
    inline virtual Bool_t IsSortable() const { return kTRUE; }
    inline Int_t GetId() const { return fId; }
    inline void SetId(Int_t Id) { fId = Id; }
    inline const TLorentzVector& P() const {
	// Return momentum 4 vector (GeV)
	return fP;
    }
    inline TVector3 P3() const {
	// Return momentum 3 vector (GeV)
	return fP.Vect();
    }
    inline Float_t Et() {
	// Return transvere energy (GeV)
	//
	// !!! The return value is cached. Do not use this function !!!
	// !!! inside TTree::Draw() or TTree::Scan().               !!!
	// !!! Use eg "fCone4H1TowerJets.fP.Et()" instead.          !!!
	//
	if ( fCompute ) ComputeTransientVars();
	return fEt;
    }
    inline Float_t Pt() {
	// Return transvere momentum (GeV)
	//
	// !!! The return value is cached. Do not use this function !!!
	// !!! inside TTree::Draw() or TTree::Scan().               !!!
	// !!! Use eg "fCone4H1TowerJets.fP.Perp()" instead.        !!!
	//
	if ( fCompute ) ComputeTransientVars();
	return fPt;
    }
    inline Float_t Eta() {
	// Return pseudo-rapidity
	//
	// !!! The return value is cached. Do not use this function !!!
	// !!! inside TTree::Draw() or TTree::Scan().               !!!
	// !!! Use eg "fCone4H1TowerJets.fP.Eta()" instead.         !!!
	//
	if ( fCompute ) ComputeTransientVars();
	return fEta;
    }
    inline Float_t Phi() {
	// Return azimuth (rad)
	//
	// !!! The return value is cached. Do not use this function !!!
	// !!! inside TTree::Draw() or TTree::Scan().               !!!
	// !!! Use eg "fCone4H1TowerJets.fP.Phi()" instead.         !!!
	//
	if ( fCompute ) ComputeTransientVars();
	return fPhi;
    }
    inline Float_t Theta() {
	// Return polar angle (rad)
	//
	// !!! The return value is cached. Do not use this function !!!
	// !!! inside TTree::Draw() or TTree::Scan().               !!!
	// !!! Use eg "fCone4H1TowerJets.fP.Theta()" instead.       !!!
	//
	if ( fCompute ) ComputeTransientVars();
	return fTheta;
    }
    inline Float_t CosTheta() {
	// Return cosine of the polar angle
	//
	// !!! The return value is cached. Do not use this function !!!
	// !!! inside TTree::Draw() or TTree::Scan().               !!!
	// !!! Use eg "fCone4H1TowerJets.fP.CosTheta()" instead.    !!!
	//
	if ( fCompute ) ComputeTransientVars();
	return fCosTheta;
    }

    inline Float_t M() const {
	// Return invariant mass (GeV)
	return fP.M();
    }
    inline Float_t Px() const {
	// Return momentum x component
	return fP.Px();
    }
    inline Float_t Py() const {
	// Return momentum y component
	return fP.Py();
    }
    inline Float_t Pz() const {
	// Return momentum z component
	return fP.Pz();
    }
    inline Float_t E() const {
	// Return energy component
	return fP.E();
    }

    inline void SetPtEtaPhiE(Float_t Pt, Float_t Eta, Float_t Phi,
			     Float_t E) {
        fP.SetPtEtaPhiE(Pt, Eta, Phi, E);
	ComputeTransientVars();
    }
    void SetP(TLorentzVector const & p) {
        fP = p;
        ComputeTransientVars();
    }
    
  private:
    void ComputeTransientVars();
    
    ClassDef(HepJet,3) // Basic jet class
};
#endif

 HepJet.h:1
 HepJet.h:2
 HepJet.h:3
 HepJet.h:4
 HepJet.h:5
 HepJet.h:6
 HepJet.h:7
 HepJet.h:8
 HepJet.h:9
 HepJet.h:10
 HepJet.h:11
 HepJet.h:12
 HepJet.h:13
 HepJet.h:14
 HepJet.h:15
 HepJet.h:16
 HepJet.h:17
 HepJet.h:18
 HepJet.h:19
 HepJet.h:20
 HepJet.h:21
 HepJet.h:22
 HepJet.h:23
 HepJet.h:24
 HepJet.h:25
 HepJet.h:26
 HepJet.h:27
 HepJet.h:28
 HepJet.h:29
 HepJet.h:30
 HepJet.h:31
 HepJet.h:32
 HepJet.h:33
 HepJet.h:34
 HepJet.h:35
 HepJet.h:36
 HepJet.h:37
 HepJet.h:38
 HepJet.h:39
 HepJet.h:40
 HepJet.h:41
 HepJet.h:42
 HepJet.h:43
 HepJet.h:44
 HepJet.h:45
 HepJet.h:46
 HepJet.h:47
 HepJet.h:48
 HepJet.h:49
 HepJet.h:50
 HepJet.h:51
 HepJet.h:52
 HepJet.h:53
 HepJet.h:54
 HepJet.h:55
 HepJet.h:56
 HepJet.h:57
 HepJet.h:58
 HepJet.h:59
 HepJet.h:60
 HepJet.h:61
 HepJet.h:62
 HepJet.h:63
 HepJet.h:64
 HepJet.h:65
 HepJet.h:66
 HepJet.h:67
 HepJet.h:68
 HepJet.h:69
 HepJet.h:70
 HepJet.h:71
 HepJet.h:72
 HepJet.h:73
 HepJet.h:74
 HepJet.h:75
 HepJet.h:76
 HepJet.h:77
 HepJet.h:78
 HepJet.h:79
 HepJet.h:80
 HepJet.h:81
 HepJet.h:82
 HepJet.h:83
 HepJet.h:84
 HepJet.h:85
 HepJet.h:86
 HepJet.h:87
 HepJet.h:88
 HepJet.h:89
 HepJet.h:90
 HepJet.h:91
 HepJet.h:92
 HepJet.h:93
 HepJet.h:94
 HepJet.h:95
 HepJet.h:96
 HepJet.h:97
 HepJet.h:98
 HepJet.h:99
 HepJet.h:100
 HepJet.h:101
 HepJet.h:102
 HepJet.h:103
 HepJet.h:104
 HepJet.h:105
 HepJet.h:106
 HepJet.h:107
 HepJet.h:108
 HepJet.h:109
 HepJet.h:110
 HepJet.h:111
 HepJet.h:112
 HepJet.h:113
 HepJet.h:114
 HepJet.h:115
 HepJet.h:116
 HepJet.h:117
 HepJet.h:118
 HepJet.h:119
 HepJet.h:120
 HepJet.h:121
 HepJet.h:122
 HepJet.h:123
 HepJet.h:124
 HepJet.h:125
 HepJet.h:126
 HepJet.h:127
 HepJet.h:128
 HepJet.h:129
 HepJet.h:130
 HepJet.h:131
 HepJet.h:132
 HepJet.h:133
 HepJet.h:134
 HepJet.h:135
 HepJet.h:136
 HepJet.h:137
 HepJet.h:138
 HepJet.h:139
 HepJet.h:140
 HepJet.h:141
 HepJet.h:142
 HepJet.h:143
 HepJet.h:144
 HepJet.h:145
 HepJet.h:146
 HepJet.h:147
 HepJet.h:148
 HepJet.h:149
 HepJet.h:150
 HepJet.h:151
 HepJet.h:152
 HepJet.h:153
 HepJet.h:154
 HepJet.h:155
 HepJet.h:156
 HepJet.h:157
 HepJet.h:158
 HepJet.h:159
 HepJet.h:160
 HepJet.h:161
 HepJet.h:162
 HepJet.h:163