ATLAS Offline Software
JetGhostThinning.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 // JetGhostThinning.cxx, (c) ATLAS Detector software
8 
11 #include "xAODBase/IParticle.h"
13 
14 // Constructor
16  const std::string &t, const std::string &n, const IInterface *p)
17  : base_class(t, n, p) {}
18 
19 // Destructor
21 
22 // Athena initialize
24  ATH_MSG_VERBOSE("initialize() ...");
25 
26  // Check that we have a jet container
27  if (m_jetSGKey.empty()) {
28  ATH_MSG_FATAL("No jet collection provided!");
29  return StatusCode::FAILURE;
30  }
31 
32  ATH_CHECK(m_jetSGKey.initialize());
33  ATH_MSG_INFO("Using " << m_jetSGKey.key() << " as the jet collection");
34 
35  // Initialize the selection string parser if provided
36  if (!m_selectionString.empty()) {
37  ATH_CHECK(initializeParser(m_selectionString));
38  ATH_MSG_INFO("Jet selection string: " << m_selectionString);
39  }
40 
41  // Check that we have a ghost name
42  if (m_ghostName.empty()) {
43  ATH_MSG_FATAL("No ghost name provided!");
44  return StatusCode::FAILURE;
45  }
46 
47  // Check that we have a ghost container name
48  if (m_ghostContainerName.empty()) {
49  ATH_MSG_FATAL("No ghost container name provided!");
50  return StatusCode::FAILURE;
51  }
52 
53  // Initialize thinning handle key
54  m_ghostContainerKey = SG::ThinningHandleKey<xAOD::IParticleContainer>(m_ghostContainerName);
55  ATH_CHECK(m_ghostContainerKey.initialize(m_streamName));
56 
57  ATH_MSG_INFO("Thinning ghost objects:");
58  ATH_MSG_INFO(" Ghost name: " << m_ghostName);
59  ATH_MSG_INFO(" Ghost container: " << m_ghostContainerName);
60 
61  return StatusCode::SUCCESS;
62 }
63 
65  ATH_MSG_VERBOSE("finalize() ...");
66  return StatusCode::SUCCESS;
67 }
68 
69 // The thinning itself
71  const EventContext &ctx = Gaudi::Hive::currentContext();
72 
73  // Retrieve jet collection
75  if (!jets.isValid()) {
76  ATH_MSG_ERROR("Failed to retrieve jet collection " << m_jetSGKey.key());
77  return StatusCode::FAILURE;
78  }
79 
80  // Get the jets that pass the selection
81  std::vector<const xAOD::Jet*> selectedJets;
82 
83  if (!m_selectionString.empty()) {
84  // Use the expression parser
85  std::vector<int> entries = m_parser->evaluateAsVector();
86  if (entries.size() != jets->size()) {
87  ATH_MSG_ERROR("Selection string evaluation size mismatch!");
88  return StatusCode::FAILURE;
89  }
90 
91  for (size_t i = 0; i < jets->size(); ++i) {
92  if (entries[i] == 1) {
93  selectedJets.push_back((*jets)[i]);
94  }
95  }
96  } else {
97  // Keep all jets
98  for (const auto* jet : *jets) {
99  selectedJets.push_back(jet);
100  }
101  }
102 
103  ATH_MSG_DEBUG("Number of selected jets: " << selectedJets.size()
104  << " out of " << jets->size());
105 
106  // Get thinning handle for the ghost container
107  SG::ThinningHandle<xAOD::IParticleContainer> ghostContainer(m_ghostContainerKey, ctx);
108 
109  if (!ghostContainer.isValid()) {
110  ATH_MSG_WARNING("Ghost container " << m_ghostContainerKey.key() << " not valid, skipping");
111  return StatusCode::SUCCESS;
112  }
113 
114  size_t nObjects = ghostContainer->size();
115  std::vector<bool> mask(nObjects, false);
116  const int maskSize = static_cast<int>(mask.size());
117 
118  // Create accessor for the ghost association
120 
121  size_t nGhostLinks = 0;
122 
123  // Loop over selected jets and mark ghost objects to keep using their indices
124  for (const auto* jet : selectedJets) {
125  if (!jet) continue;
126 
127  // Check if this jet has the ghost association
128  if (!ghostAcc.isAvailable(*jet)) {
129  ATH_MSG_VERBOSE("Jet does not have ghost association: " << m_ghostName);
130  continue;
131  }
132 
133  // Get the ghost links
134  const std::vector<ElementLink<xAOD::IParticleContainer>>& ghostLinks = ghostAcc(*jet);
135 
136  ATH_MSG_VERBOSE("Jet has " << ghostLinks.size() << " ghost objects");
137 
138  // Mark all valid ghost objects using their indices
139  for (const auto& link : ghostLinks) {
140  if (!link.isValid()) continue;
141 
142  int index = link.index();
143  if (index < 0) {
144  ATH_MSG_VERBOSE(" Ghost link has negative index, skipping");
145  continue;
146  }
147 
148  if (index < maskSize) {
149  mask[index] = true;
150  nGhostLinks++;
151  const xAOD::IParticle* ghostObj = *link;
152  if (ghostObj) {
153  ATH_MSG_VERBOSE(" Ghost object at index " << index << ": pt=" << ghostObj->pt()
154  << ", eta=" << ghostObj->eta());
155  }
156  }
157  }
158  }
159 
160  size_t nKept = std::count(mask.begin(), mask.end(), true);
161  ATH_MSG_DEBUG("Found " << nGhostLinks << " ghost links from selected jets");
162  ATH_MSG_DEBUG("Ghost container " << m_ghostContainerKey.key() << ": keeping "
163  << nKept << " out of " << nObjects << " objects");
164 
165  ghostContainer.keep(mask);
166 
167  return StatusCode::SUCCESS;
168 }
ATH_MSG_FATAL
#define ATH_MSG_FATAL(x)
Definition: AthMsgStreamMacros.h:34
IParticle.h
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
DerivationFramework::JetGhostThinning::doThinning
virtual StatusCode doThinning() const override
Definition: JetGhostThinning.cxx:70
SG::ThinningHandleKey
HandleKey object for adding thinning to an object.
Definition: ThinningHandleKey.h:38
ThinningHandle.h
Handle for requesting thinning for a data object.
SG::ReadHandle
Definition: StoreGate/StoreGate/ReadHandle.h:67
index
Definition: index.py:1
defineDB.jets
jets
Definition: JetTagCalibration/share/defineDB.py:24
SG::ConstAccessor
Helper class to provide constant type-safe access to aux data.
Definition: ConstAccessor.h:55
read_hist_ntuple.t
t
Definition: read_hist_ntuple.py:5
ATH_MSG_VERBOSE
#define ATH_MSG_VERBOSE(x)
Definition: AthMsgStreamMacros.h:28
xAOD::IParticle
Class providing the definition of the 4-vector interface.
Definition: Event/xAOD/xAODBase/xAODBase/IParticle.h:41
XMLtoHeader.count
count
Definition: XMLtoHeader.py:84
python.utils.AtlRunQueryLookup.mask
string mask
Definition: AtlRunQueryLookup.py:459
SG::ThinningHandle
Handle for requesting thinning for a data object.
Definition: ThinningHandle.h:84
SG::ThinningHandleBase::keep
void keep(size_t ndx)
Mark that index ndx in the container should be kept (not thinned away).
Definition: ThinningHandleBase.cxx:75
JetGhostThinning.h
ROBbits::maskSize
const uint32_t maskSize
Definition: TrigMonROBData.cxx:17
python.utils.AtlRunQueryDQUtils.p
p
Definition: AtlRunQueryDQUtils.py:209
jet
Definition: JetCalibTools_PlotJESFactors.cxx:23
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
lumiFormat.i
int i
Definition: lumiFormat.py:85
beamspotman.n
n
Definition: beamspotman.py:727
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
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
DerivationFramework::JetGhostThinning::~JetGhostThinning
virtual ~JetGhostThinning()
DerivationFramework::JetGhostThinning::JetGhostThinning
JetGhostThinning(const std::string &t, const std::string &n, const IInterface *p)
Definition: JetGhostThinning.cxx:15
xAOD::IParticle::pt
virtual double pt() const =0
The transverse momentum ( ) of the particle.
DerivationFramework::JetGhostThinning::initialize
virtual StatusCode initialize() override
Definition: JetGhostThinning.cxx:23
DeMoScan.index
string index
Definition: DeMoScan.py:362
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
xAOD::IParticle::eta
virtual double eta() const =0
The pseudorapidity ( ) of the particle.
entries
double entries
Definition: listroot.cxx:49
SG::ConstAccessor::isAvailable
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
DerivationFramework::JetGhostThinning::finalize
virtual StatusCode finalize() override
Definition: JetGhostThinning.cxx:64
AuxElement.h
Base class for elements of a container that can have aux data.