ATLAS Offline Software
JetReclusteringAlg.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 
8 #include "Math/Vector4D.h"
9 #include "Math/VectorUtil.h"
10 
11 namespace CP {
12  using ROOT::Math::PtEtaPhiEVector;
14 
16 
19 
22 
24 
26 
27  // reclustering properties
28  if (m_reclusteredJetsRadius < 0) {
29  ANA_MSG_ERROR("reclusteredJetsRadius is negative, please fix");
30  return StatusCode::FAILURE;
31  }
32 
33  fastjet::JetAlgorithm algorithm;
34  static const std::unordered_map<std::string, fastjet::JetAlgorithm> algorithmMap = {
35  {"AntiKt", fastjet::antikt_algorithm},
36  {"Kt", fastjet::kt_algorithm},
38  };
39  auto it = algorithmMap.find(m_clusteringAlgorithm.value());
40  if (it != algorithmMap.end()) {
41  algorithm = it->second;
42  } else {
43  ATH_MSG_ERROR("Jet reclustering algorithm " << m_clusteringAlgorithm.value()
44  << " is not supported by CP::JetReclusteringAlg! Choose amongst: AntiKt, Kt, CamKt.");
45  return StatusCode::FAILURE;
46  }
47 
48  m_fastjetClustering = std::make_unique<fastjet::JetDefinition>(algorithm, m_reclusteredJetsRadius, fastjet::E_scheme, fastjet::Best);
49 
50  return StatusCode::SUCCESS;
51  }
52 
54 
55  for (const auto &sys : m_systematicsList.systematicsVector()) {
56 
57  const xAOD::JetContainer *jets = nullptr;
59 
60  auto rcJetAuxContainer = std::make_unique<xAOD::JetAuxContainer>();
61  auto rcJetContainer = std::make_unique<xAOD::JetContainer>();
62  rcJetContainer->setStore(rcJetAuxContainer.get());
63 
64  // only take selected jets
65  std::vector<fastjet::PseudoJet> clusters;
66  std::vector<const xAOD::Jet*> selectedJets;
67  for (const xAOD::Jet *jet : *jets) {
68  if (m_jetSelection.getBool(*jet, sys)) {
69  fastjet::PseudoJet pseudoJet(jet->px(), jet->py(), jet->pz(), jet->e());
70  clusters.emplace_back(std::move(pseudoJet));
71  selectedJets.push_back(jet);
72  }
73  }
74 
75  // no jets (clusters) - nothing to do
76  if (clusters.empty()) {
77  ANA_CHECK(m_outHandle.record(std::move(rcJetContainer), std::move(rcJetAuxContainer), sys));
78  continue;
79  }
80 
82  std::vector<fastjet::PseudoJet> rcJets = cluster_sequence.inclusive_jets(0.0);
83 
84  for (const auto& irc : rcJets) {
85  xAOD::Jet *rcJet = new xAOD::Jet();
86 
87  const xAOD::JetFourMom_t p4(irc.pt(), irc.eta(), irc.phi_std(), irc.m());
88  rcJetContainer->push_back(rcJet);
89 
90  rcJetContainer->back()->setJetP4(p4);
91 
92  const auto constituents = irc.constituents();
93 
94  const std::vector<int> indices = this->matchRCjets(selectedJets, constituents);
95 
96  if (indices.size() != constituents.size()) {
97  ANA_MSG_ERROR("Size of the contituents does not match the size of the indices");
98  ANA_CHECK(m_outHandle.record(std::move(rcJetContainer), std::move(rcJetAuxContainer), sys));
99  return StatusCode::FAILURE;
100  }
101 
102  m_smallRjetIndicesDecor.set(*(rcJetContainer->back()), indices, sys);
103 
104  // Write energy to the RC jet, since default 4-vector form is (pt, eta, phi, m)
105  m_rcEnergyDecor.set(*(rcJetContainer->back()), irc.e(), sys);
106  }
107 
108  ANA_CHECK(m_outHandle.record(std::move(rcJetContainer), std::move(rcJetAuxContainer), sys));
109 
110  }
111  return StatusCode::SUCCESS;
112  }
113 
114  std::vector<int> JetReclusteringAlg::matchRCjets(const std::vector<const xAOD::Jet*>& smallJets, const std::vector<fastjet::PseudoJet>& constituents) const {
115  std::vector<int> indices;
116  for (const auto& iconst : constituents) {
117  PtEtaPhiEVector iconstP4(iconst.pt(), iconst.eta(), iconst.phi_std(), iconst.e());
118  for (std::size_t i = 0; i < smallJets.size(); ++i) {
119  const float dR = DeltaR(smallJets.at(i)->p4(), iconstP4);
120  if (dR > 0.1) continue; // no match
121 
122  indices.emplace_back(i);
123  }
124  }
125  return indices;
126  }
127 
128 } // namespace
CP::SysWriteDecorHandle::set
void set(const SG::AuxElement &object, const T &value, const CP::SystematicSet &sys) const
set the object decoration for the given systematic
algorithm
std::string algorithm
Definition: hcg.cxx:82
CP::JetReclusteringAlg::initialize
virtual StatusCode initialize() override
Definition: JetReclusteringAlg.cxx:15
xAOD::JetAlgorithmType::kt_algorithm
@ kt_algorithm
Definition: JetContainerInfo.h:31
CP::JetReclusteringAlg::m_reclusteredJetsRadius
Gaudi::Property< float > m_reclusteredJetsRadius
Definition: JetReclusteringAlg.h:62
met::DeltaR
@ DeltaR
Definition: METRecoCommon.h:11
Trk::indices
std::pair< long int, long int > indices
Definition: AlSymMatBase.h:24
skel.it
it
Definition: skel.GENtoEVGEN.py:396
ANA_MSG_ERROR
#define ANA_MSG_ERROR(xmsg)
Macro printing error messages.
Definition: Control/AthToolSupport/AsgMessaging/AsgMessaging/MessageCheck.h:294
ANA_CHECK
#define ANA_CHECK(EXP)
check whether the given expression was successful
Definition: Control/AthToolSupport/AsgMessaging/AsgMessaging/MessageCheck.h:324
CP::SysReadHandle::retrieve
::StatusCode retrieve(const T *&object, const CP::SystematicSet &sys) const
retrieve the object for the given name
CP::SysListHandle::systematicsVector
const std::vector< CP::SystematicSet > & systematicsVector() const
the list of systematics to loop over
Definition: SysListHandle.cxx:96
CP
Select isolated Photons, Electrons and Muons.
Definition: Control/xAODRootAccess/xAODRootAccess/TEvent.h:49
mapkey::sys
@ sys
Definition: TElectronEfficiencyCorrectionTool.cxx:42
CP::JetReclusteringAlg::execute
virtual StatusCode execute() override
Definition: JetReclusteringAlg.cxx:53
JetReclusteringAlg.h
CP::SysReadHandle::initialize
StatusCode initialize(SysListHandle &sysListHandle)
initialize this handle
jet
Definition: JetCalibTools_PlotJESFactors.cxx:23
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
CP::SysListHandle::initialize
::StatusCode initialize()
intialize this property
Definition: SysListHandle.cxx:69
CP::SysReadSelectionHandle::getBool
bool getBool(const SG::AuxElement &element, const CP::SystematicSet &sys) const
get the selection as a bool
jet::ClusterSequence
fastjet::ClusterSequence ClusterSequence
Definition: ClusterSequence.h:21
lumiFormat.i
int i
Definition: lumiFormat.py:85
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
CP::SysWriteDecorHandle::initialize
StatusCode initialize(SysListHandle &sysListHandle, const ISysHandleBase &objectHandle)
initialize this handle
CP::JetReclusteringAlg::m_jetsHandle
CP::SysReadHandle< xAOD::JetContainer > m_jetsHandle
Definition: JetReclusteringAlg.h:38
CP::JetReclusteringAlg::m_clusteringAlgorithm
Gaudi::Property< std::string > m_clusteringAlgorithm
Definition: JetReclusteringAlg.h:59
CP::JetReclusteringAlg::m_rcEnergyDecor
CP::SysWriteDecorHandle< float > m_rcEnergyDecor
Definition: JetReclusteringAlg.h:49
CP::JetReclusteringAlg::m_fastjetClustering
std::unique_ptr< fastjet::JetDefinition > m_fastjetClustering
Definition: JetReclusteringAlg.h:65
DataVector
Derived DataVector<T>.
Definition: DataVector.h:794
xAOD::JetFourMom_t
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< double > > JetFourMom_t
Base 4 Momentum type for Jet.
Definition: JetTypes.h:17
CP::JetReclusteringAlg::m_systematicsList
CP::SysListHandle m_systematicsList
Definition: JetReclusteringAlg.h:36
CP::JetReclusteringAlg::m_jetSelection
CP::SysReadSelectionHandle m_jetSelection
Definition: JetReclusteringAlg.h:41
xAOD::JetAlgorithmType::cambridge_algorithm
@ cambridge_algorithm
Definition: JetContainerInfo.h:32
xAOD::JetAlgorithmType::antikt_algorithm
@ antikt_algorithm
Definition: JetContainerInfo.h:33
CP::JetReclusteringAlg::matchRCjets
std::vector< int > matchRCjets(const std::vector< const xAOD::Jet * > &smallJets, const std::vector< fastjet::PseudoJet > &constituents) const
Definition: JetReclusteringAlg.cxx:114
xAOD::Jet_v1
Class describing a jet.
Definition: Jet_v1.h:57
CP::SysReadSelectionHandle::initialize
StatusCode initialize(SysListHandle &sysListHandle, const ISysHandleBase &objectHandle)
initialize the accessor
Definition: SysReadSelectionHandle.cxx:34
RunTileMonitoring.clusters
clusters
Definition: RunTileMonitoring.py:133
defineDB.jets
list jets
Definition: JetTagCalibration/share/defineDB.py:24
CP::JetReclusteringAlg::m_outHandle
CP::SysWriteHandle< xAOD::JetContainer, xAOD::JetAuxContainer > m_outHandle
Definition: JetReclusteringAlg.h:54
SG::AllowEmpty
@ AllowEmpty
Definition: StoreGate/StoreGate/VarHandleKey.h:30
CP::JetReclusteringAlg::m_smallRjetIndicesDecor
CP::SysWriteDecorHandle< std::vector< int > > m_smallRjetIndicesDecor
Definition: JetReclusteringAlg.h:45
xAOD::Jet
Jet_v1 Jet
Definition of the current "jet version".
Definition: Event/xAOD/xAODJet/xAODJet/Jet.h:17