ATLAS Offline Software
FTAGValidationAlgorithm.cxx
Go to the documentation of this file.
1 //Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
2 
4 
5 namespace FTAGValidation {
6 
8  ISvcLocator* pSvcLocator ) :
9  ::AthAlgorithm( name, pSvcLocator ) {}
10 
12 
14  ATH_MSG_INFO( "Inizializing " << name() << " ..." );
16  return StatusCode::SUCCESS;
17  }
18 
19  /* ===================================================================== */
20 
22 
23  if ( jet->p4().Et() / Gaudi::Units::GeV < m_minJetPt ) return false;
24  if ( std::fabs( jet->eta() ) < m_minJetEta ) return false;
25  if ( std::fabs( jet->eta() ) > m_maxJetEta ) return false;
26 
27  return true;
28  }
29 
31 
32  if ( jet->p4().Et() / Gaudi::Units::GeV < m_minTrigJetPt ) return false;
33  if ( std::fabs( jet->eta() ) > m_maxTrigJetEta ) return false;
34 
35  return true;
36  }
37 
39 
40  // This part is inspired from 'PhysicsAnalysis/JetTagging/JetTagMonitoring/src/JetTagMonitorAlgorithm.cxx'
41 
42  float jetQuality = 0;
43  if ( !jet->getAttribute(xAOD::JetAttribute::LArQuality,jetQuality) ) {
44  ATH_MSG_ERROR( "JetAttribute 'xAOD::JetAttribute::LArQuality' is not available." );
45  return false;
46  }
47 
48  float jetTime = 0;
49  if ( !jet->getAttribute(xAOD::JetAttribute::Timing,jetTime) ) {
50  ATH_MSG_ERROR( "JetAttribute 'xAOD::JetAttribute::Timing' is not available." );
51  return false;
52  }
53 
54  float hecq = 0;
55  if ( !jet->getAttribute(xAOD::JetAttribute::HECQuality,hecq) ) {
56  ATH_MSG_ERROR( "JetAttribute 'xAOD::JetAttribute::HECQuality' is not available." );
57  return false;
58  }
59 
60  float negE = 0;
61  if ( !jet->getAttribute(xAOD::JetAttribute::NegativeE,negE) ) {
62  ATH_MSG_ERROR( "JetAttribute 'xAOD::JetAttribute::NegativeE' is not available." );
63  return false;
64  }
65 
66  std::vector< float > SumPtTrkPt1000;
67  if ( !jet->getAttribute(xAOD::JetAttribute::SumPtTrkPt1000,SumPtTrkPt1000) ) {
68  ATH_MSG_ERROR( "JetAttribute 'xAOD::JetAttribute::SumPtTrkPt1000' is not available." );
69  return false;
70  }
71  float chf = SumPtTrkPt1000.size() > 0 ? SumPtTrkPt1000.at(0)/jet->pt() : -1;
72 
73  float emf = 0;
74  if ( !jet->getAttribute(xAOD::JetAttribute::EMFrac,emf) ) {
75  ATH_MSG_ERROR( "JetAttribute 'xAOD::JetAttribute::EMFrac' is not available." );
76  return false;
77  }
78 
79  float hecf = 0;
80  if ( !jet->getAttribute(xAOD::JetAttribute::HECFrac,hecf) ) {
81  ATH_MSG_ERROR( "JetAttribute 'xAOD::JetAttribute::HECFrac' is not available." );
82  return false;
83  }
84 
85  float fracSamplingMax = 0;
86  if ( !jet->getAttribute(xAOD::JetAttribute::FracSamplingMax,fracSamplingMax) ) {
87  ATH_MSG_ERROR( "JetAttribute 'xAOD::JetAttribute::FracSamplingMax' is not available." );
88  return false;
89  }
90 
91  if ( hecf > 0.5 && fabs( hecq ) > 0.5 ) return false;
92  if ( fabs( negE ) > 60*Gaudi::Units::GeV ) return false;
93  if ( emf > 0.95 && fabs( jetQuality ) > 0.8 && fabs( jet->eta() ) < 2.8 ) return false;
94  if ( fabs( jetTime ) > 25 ) return false;
95  if ( emf < 0.05 && chf < 0.05 && fabs( jet->eta() ) < 2 ) return false;
96  if ( emf< 0.05 && fabs( jet->eta() ) >= 2 ) return false;
97  if ( fracSamplingMax > 0.99 && fabs( jet->eta() ) < 2 ) return false;
98 
99  return true;
100  }
101 
103  const std::string& jetType ) const {
104 
105  if ( jet->isAvailable< float >( "Jvt" ) == false ) {
106  ATH_MSG_ERROR( "Jet auxdata 'Jvt' is not available." );
107  return false;
108  }
109 
110  float jvt = jet->auxdata< float >( "Jvt" );
111  double jvtCut = 0.59;
112  if ( jetType == "AntiKt4EMPFlowJets" )
113  jvtCut = 0.2;
114 
115  if ( jet->pt() / Gaudi::Units::GeV >= 60. )
116  return true;
117  if ( fabs( jet->eta() ) >= 2.4 )
118  return true;
119  if ( jvt > jvtCut )
120  return true;
121 
122  return false;
123  }
124 
126 
127  if ( vertexContainer->size() == 0 ) {
128  ATH_MSG_WARNING( "Vertex Container has size 0! This can't be right!" );
129  return nullptr;
130  }
131 
132  for ( unsigned int i(0); i<vertexContainer->size(); i++ ) {
133  const xAOD::Vertex *vertex = vertexContainer->at(i);
134  if ( vertex->vertexType() != xAOD::VxType::VertexType::PriVtx ) continue;
135  return vertex;
136  }
137 
138  ATH_MSG_DEBUG( "None of the vertexes in the vertex container is a primary vertex!" );
139  ATH_MSG_DEBUG( "Using Dummy Vertex." );
140  return vertexContainer->at(0);
141  }
142 
143  int FTAGValidationAlgorithm::getMatchedOfflineJetIndex( const xAOD::Jet* onJet, std::vector< const xAOD::Jet* > offJets ) const {
144  int matchedOffJetIndex = -1;
145  double minDr = 0;
146 
147  TLorentzVector onJet_4vec, offJet_4vec;
148  onJet_4vec.SetPtEtaPhiM(onJet->pt(), onJet->eta(), onJet->phi(), onJet->m());
149 
150  for ( unsigned int index=0; index<offJets.size(); index++ ) {
151  const xAOD::Jet *offJet = offJets.at( index );
152 
153  offJet_4vec.SetPtEtaPhiM(offJet->pt(), offJet->eta(), offJet->phi(), offJet->m());
154  double deltaR = onJet_4vec.DeltaR( offJet_4vec );
155  ATH_MSG_DEBUG("deltaR: " << deltaR );
156 
157  if ( deltaR > 0.1 ) continue;
158  if ( deltaR < minDr || matchedOffJetIndex == -1 ) {
159  minDr = deltaR;
160  matchedOffJetIndex = index;
161  }
162  }
163  ATH_MSG_DEBUG("minDr: " << minDr );
164  ATH_MSG_DEBUG("matchedOffJetIndex: " << matchedOffJetIndex );
165  return matchedOffJetIndex;
166  }
167 
168 }
169 
170 
171 
172 
FTAGValidation::FTAGValidationAlgorithm::passTriggerJetKinematicSelection
bool passTriggerJetKinematicSelection(const xAOD::Jet *) const
Definition: FTAGValidationAlgorithm.cxx:30
FTAGValidation::FTAGValidationAlgorithm::m_minJetEta
Gaudi::Property< float > m_minJetEta
Definition: FTAGValidationAlgorithm.h:50
xAOD::JetAttribute::LArQuality
@ LArQuality
Definition: JetAttributes.h:83
FTAGValidation::FTAGValidationAlgorithm::m_minTrigJetPt
Gaudi::Property< float > m_minTrigJetPt
Definition: FTAGValidationAlgorithm.h:62
xAOD::JetAttribute::FracSamplingMax
@ FracSamplingMax
Definition: JetAttributes.h:116
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
FTAGValidation::FTAGValidationAlgorithm::FTAGValidationAlgorithm
FTAGValidationAlgorithm()
index
Definition: index.py:1
FTAGValidation::FTAGValidationAlgorithm::passJetQualitySelection
bool passJetQualitySelection(const xAOD::Jet *) const
Definition: FTAGValidationAlgorithm.cxx:38
xAOD::Jet_v1::phi
virtual double phi() const
The azimuthal angle ( ) of the particle.
Definition: Jet_v1.cxx:54
xAOD::JetAttribute::HECQuality
@ HECQuality
Definition: JetAttributes.h:80
FTAGValidationAlgorithm.h
jet
Definition: JetCalibTools_PlotJESFactors.cxx:23
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
FTAGValidation
Definition: FTAGValidationAlgorithm.h:16
lumiFormat.i
int i
Definition: lumiFormat.py:92
xAOD::JetAttribute::EMFrac
@ EMFrac
Definition: JetAttributes.h:112
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
xAOD::JetAttribute::SumPtTrkPt1000
@ SumPtTrkPt1000
Definition: JetAttributes.h:107
xAOD::VxType::PriVtx
@ PriVtx
Primary vertex.
Definition: TrackingPrimitives.h:571
FTAGValidation::FTAGValidationAlgorithm::m_maxTrigJetEta
Gaudi::Property< float > m_maxTrigJetEta
Definition: FTAGValidationAlgorithm.h:59
xAOD::JetAttribute::Timing
@ Timing
Definition: JetAttributes.h:90
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
FTAGValidation::FTAGValidationAlgorithm::m_minJetPt
Gaudi::Property< float > m_minJetPt
Definition: FTAGValidationAlgorithm.h:56
SG::VarHandleKey::initialize
StatusCode initialize(bool used=true)
If this object is used as a property, then this should be called during the initialize phase.
Definition: AthToolSupport/AsgDataHandles/Root/VarHandleKey.cxx:103
FTAGValidation::FTAGValidationAlgorithm::m_eventInfoKey
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfoKey
Definition: FTAGValidationAlgorithm.h:47
DataVector
Derived DataVector<T>.
Definition: DataVector.h:581
AthAlgorithm
Definition: AthAlgorithm.h:47
xAOD::Jet_v1::eta
virtual double eta() const
The pseudorapidity ( ) of the particle.
Definition: Jet_v1.cxx:49
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
FTAGValidation::FTAGValidationAlgorithm::passJetJVTSelection
bool passJetJVTSelection(const xAOD::Jet *, const std::string &jetType="AntiKt4EMTopoJets") const
Definition: FTAGValidationAlgorithm.cxx:102
FTAGValidation::FTAGValidationAlgorithm::~FTAGValidationAlgorithm
virtual ~FTAGValidationAlgorithm()
Definition: FTAGValidationAlgorithm.cxx:11
xAOD::Jet_v1
Class describing a jet.
Definition: Jet_v1.h:57
Trk::vertex
@ vertex
Definition: MeasurementType.h:21
DeMoScan.index
string index
Definition: DeMoScan.py:362
ConvertOldUJHistosToNewHistos.jetType
string jetType
Definition: ConvertOldUJHistosToNewHistos.py:121
xAOD::Jet_v1::m
virtual double m() const
The invariant mass of the particle.
Definition: Jet_v1.cxx:59
xAOD::Vertex_v1
Class describing a Vertex.
Definition: Vertex_v1.h:42
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
xAOD::JetAttribute::HECFrac
@ HECFrac
Definition: JetAttributes.h:113
DataVector::at
const T * at(size_type n) const
Access an element, as an rvalue.
FTAGValidation::FTAGValidationAlgorithm::m_maxJetEta
Gaudi::Property< float > m_maxJetEta
Definition: FTAGValidationAlgorithm.h:53
xAOD::Jet_v1::pt
virtual double pt() const
The transverse momentum ( ) of the particle.
Definition: Jet_v1.cxx:44
GeV
#define GeV
Definition: CaloTransverseBalanceVecMon.cxx:30
FTAGValidation::FTAGValidationAlgorithm::passJetKinematicSelection
bool passJetKinematicSelection(const xAOD::Jet *) const
Definition: FTAGValidationAlgorithm.cxx:21
makeComparison.deltaR
float deltaR
Definition: makeComparison.py:36
DataVector::size
size_type size() const noexcept
Returns the number of elements in the collection.
FTAGValidation::FTAGValidationAlgorithm::initialize
virtual StatusCode initialize()
Definition: FTAGValidationAlgorithm.cxx:13
FTAGValidation::FTAGValidationAlgorithm::getPrimaryVertex
const xAOD::Vertex * getPrimaryVertex(const xAOD::VertexContainer *) const
Definition: FTAGValidationAlgorithm.cxx:125
xAOD::JetAttribute::NegativeE
@ NegativeE
Definition: JetAttributes.h:84
FTAGValidation::FTAGValidationAlgorithm::getMatchedOfflineJetIndex
int getMatchedOfflineJetIndex(const xAOD::Jet *, std::vector< const xAOD::Jet * >) const
Definition: FTAGValidationAlgorithm.cxx:143