ATLAS Offline Software
Loading...
Searching...
No Matches
HIJetRecDefs.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef HIJETREC_HIJETRECDEFS_H
6#define HIJETREC_HIJETRECDEFS_H
7
9#include "xAODJet/JetTypes.h"
12#include "TLorentzVector.h"
13
14
15namespace HIJetRec{
16
17 //define conventions for HIJets in terms of various signal states
18 constexpr const char* unsubtractedJetState() { return "JetUnsubtractedScaleMomentum";}
19 constexpr const char* subtractedJetState() { return "JetSubtractedScaleMomentum";}
20 constexpr const char* subtractedOriginCorrectedJetState() { return "JetSubtractedOriginCorrectedScaleMomentum";}
21
24
28
29 inline bool inTowerBoundary(float eta0, float phi0, float eta, float phi)
30 {
31 if( 2.*std::abs(eta-eta0) > HI::TowerBins::getBinSizeEta() ) return false;
32 if( 2.*std::abs(xAOD::P4Helpers::deltaPhi(phi0,phi)) > HI::TowerBins::getBinSizePhi() ) return false;
33 return true;
34 }
35
37 {
38 float energy=p.Energy();
39 float eta=p.Eta();
40 float phi=p.Phi();
41 if(energy < 0.)
42 {
43 eta*=-1.;
44 if(phi > 0.) phi-=M_PI;
45 else phi+=M_PI;
46 }
48 {
49 cl->setRawE(energy);
50 cl->setRawEta(eta);
51 cl->setRawPhi(phi);
52 cl->setRawM(0);
53 }
55 {
56 cl->setAltE(energy);
57 cl->setAltEta(eta);
58 cl->setAltPhi(phi);
59 cl->setAltM(0);
60 }
62 {
63 cl->setCalE(energy);
64 cl->setCalEta(eta);
65 cl->setCalPhi(phi);
66 cl->setCalM(0);
67 }
68 }
69
73 inline TLorentzVector getClusterP4( const xAOD::CaloCluster* cl, xAOD::CaloCluster::State s ){
74 TLorentzVector correct_clusP4;
75 float energy=cl->e(s);
76 float eta=cl->eta(s);
77 float phi=cl->phi(s);
78 if(energy < 0.)
79 {
80 eta*=-1.;
81 if(phi > 0.) phi-=M_PI;
82 else phi+=M_PI;
83 }
84 float pt=energy/std::cosh(eta);
85 correct_clusP4.SetPtEtaPhiE(pt, eta,phi ,energy);
86 return correct_clusP4;
87 }
88}
89#endif
#define M_PI
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
State
enum of possible signal states.
IParticle::FourMom_t FourMom_t
Definition of the 4-momentum type.
constexpr const char * unsubtractedJetState()
constexpr xAOD::JetConstitScale subtractedOriginCorrectedConstitState()
constexpr const char * subtractedJetState()
TLorentzVector getClusterP4(const xAOD::CaloCluster *cl, xAOD::CaloCluster::State s)
: getClusterP4 Added during Rel22 migration.
bool inTowerBoundary(float eta0, float phi0, float eta, float phi)
constexpr const char * subtractedOriginCorrectedJetState()
constexpr xAOD::CaloCluster::State unsubtractedClusterState()
constexpr xAOD::CaloCluster::State subtractedClusterState()
constexpr xAOD::JetConstitScale subtractedConstitState()
void setClusterP4(const xAOD::CaloCluster::FourMom_t &p, xAOD::CaloCluster *cl, xAOD::CaloCluster::State s)
constexpr xAOD::CaloCluster::State subtractedOriginCorrectedClusterState()
constexpr float getBinSizeEta()
Definition HIEventDefs.h:31
constexpr float getBinSizePhi()
Definition HIEventDefs.h:32
double deltaPhi(double phiA, double phiB)
delta Phi in range [-pi,pi[
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
JetConstitScale
Definition JetTypes.h:20
@ UncalibratedJetConstituent
Definition JetTypes.h:21
@ CalibratedJetConstituent
Definition JetTypes.h:22