ATLAS Offline Software
JetTrackParticleThinning.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 // JetTrackParticleThinning.cxx, (c) ATLAS Detector software
8 // Author: James Catmore (James.Catmore@cern.ch)
9 // This is a trivial example of an implementation of a thinning tool
10 // which removes all ID tracks which do not pass a user-defined cut
11 
13 #include "xAODJet/JetContainer.h"
16 #include "GaudiKernel/ThreadLocalContext.h"
17 #include <vector>
18 #include <string>
19 
20 // Constructor
22  const std::string& n,
23  const IInterface* p ) :
24 base_class(t,n,p)
25 {
26 }
27 
28 // Destructor
30 
31 // Athena initialize and finalize
33 {
34  // Decide which collections need to be checked for ID TrackParticles
35  ATH_MSG_VERBOSE("initialize() ...");
36  ATH_CHECK( m_inDetSGKey.initialize (m_streamName) );
37  ATH_MSG_INFO("Using " << m_inDetSGKey << "as the source collection for inner detector track particles");
38  if (m_jetKey.key().empty()) {
39  ATH_MSG_FATAL("No jet collection provided for thinning.");
40  return StatusCode::FAILURE;
41  }
42  ATH_CHECK( m_jetKey.initialize() );
43  ATH_MSG_INFO("Inner detector track particles associated with objects in " << m_jetKey.key() << " will be retained in this format with the rest being thinned away");
44 
45  // Set up the text-parsing machinery for selectiong the jet directly according to user cuts
46  if (!m_selectionString.empty() || !m_trackSelectionString.empty()) {
47  // order must match enum order EJetTrPThinningParser
48  ATH_CHECK( initializeParser( {m_selectionString, m_trackSelectionString } ));
49  }
50  return StatusCode::SUCCESS;
51 }
52 
54 {
55  ATH_MSG_VERBOSE("finalize() ...");
56  ATH_MSG_INFO("Processed "<< m_ntot <<" tracks, "<< m_npass<< " were retained ");
57  return StatusCode::SUCCESS;
58 }
59 
60 // The thinning itself
62 {
63  const EventContext& ctx = Gaudi::Hive::currentContext();
64 
65  // Retrieve main TrackParticle collection
67  (m_inDetSGKey, ctx);
68 
69  // Check the event contains tracks
70  unsigned int nTracks = importedTrackParticles->size();
71  if (nTracks==0) return StatusCode::SUCCESS;
72 
73  // Set up a mask with the same entries as the full TrackParticle collection
74  std::vector<bool> mask;
75  mask.assign(nTracks,false); // default: don't keep any tracks
76  m_ntot += nTracks;
77 
78  // Retrieve containers
79  // ... jets
80  SG::ReadHandle<xAOD::JetContainer> importedJets(m_jetKey,ctx);
81  if (!importedJets.isValid()) {
82  ATH_MSG_ERROR("No jet collection with name " << m_jetKey.key() << " found in StoreGate!");
83  return StatusCode::FAILURE;
84  }
85  unsigned int nJets(importedJets->size());
86  std::vector<const xAOD::Jet*> jetToCheck; jetToCheck.clear();
87 
88  // Execute the text parser if requested
89  if (!m_selectionString.empty()) {
90  std::vector<int> entries = m_parser[kJetSelection]->evaluateAsVector();
91  unsigned int nEntries = entries.size();
92  // check the sizes are compatible
93  if (nJets != nEntries ) {
94  ATH_MSG_ERROR("Sizes incompatible! Are you sure your selection string used jets??");
95  return StatusCode::FAILURE;
96  } else {
97  // identify which jets to keep for the thinning check
98  for (unsigned int i=0; i<nJets; ++i) if (entries[i]==1) jetToCheck.push_back((*importedJets)[i]);
99  }
100  }
101 
102  // Set elements in the mask to true if there is a corresponding ElementLink from a reconstructed object
103  // ... jets
104  if (m_selectionString.empty()) { // check all jets as user didn't provide a selection string
105  for (const auto *jetIt : *importedJets) {
106  std::vector<const xAOD::TrackParticle*> jetTracks;
107  bool haveJetTracks = jetIt->getAssociatedObjects(xAOD::JetAttribute::GhostTrack, jetTracks);
108  if ( !haveJetTracks ) {ATH_MSG_WARNING("Associated tracks not found");}
109  else {
110  for (auto & jetTrack : jetTracks) {
111  int index = jetTrack->index();
112  mask[index] = true;
113  }
114  }
115  }
116  } else { // check only jets passing user selection string
117  for (auto & jetIt : jetToCheck) {
118  std::vector<const xAOD::TrackParticle*> jetTracks;
119  bool haveJetTracks = jetIt->getAssociatedObjects(xAOD::JetAttribute::GhostTrack, jetTracks);
120  if ( !haveJetTracks ) {ATH_MSG_WARNING("Associated tracks not found");}
121  else {
122  for (auto & jetTrack : jetTracks) {
123  int index = jetTrack->index();
124  mask[index] = true;
125  }
126  }
127  }
128  }
129 
130  // Apply a track selection string.
131  if (!m_trackSelectionString.empty()) {
132  std::vector<int> entries = m_parser[kTrackThinning]->evaluateAsVector();
133  unsigned int nEntries = entries.size();
134  // check the sizes are compatible
135  if (nTracks != nEntries ) {
136  ATH_MSG_ERROR("Sizes incompatible! Are you sure your track selection string used tracks??");
137  return StatusCode::FAILURE;
138  } else {
139  // identify which jets to keep for the thinning check
140  for (unsigned int i=0; i<nEntries; ++i) {
141  if (!entries[i]) mask[i] = false;
142  }
143  }
144  }
145 
146  // Count up the mask contents
147  unsigned int n_pass=0;
148  for (unsigned int i=0; i<nTracks; ++i) {
149  if (mask[i]) ++n_pass;
150  }
151  m_npass += n_pass;
152 
153  // Execute the thinning service based on the mask. Finish.
154  importedTrackParticles.keep (mask);
155 
156  return StatusCode::SUCCESS;
157 }
158 
ATH_MSG_FATAL
#define ATH_MSG_FATAL(x)
Definition: AthMsgStreamMacros.h:34
python.PerfMonSerializer.p
def p
Definition: PerfMonSerializer.py:743
ATH_MSG_INFO
#define ATH_MSG_INFO(x)
Definition: AthMsgStreamMacros.h:31
ThinningHandle.h
Handle for requesting thinning for a data object.
SG::ReadHandle
Definition: StoreGate/StoreGate/ReadHandle.h:70
index
Definition: index.py:1
read_hist_ntuple.t
t
Definition: read_hist_ntuple.py:5
ATH_MSG_VERBOSE
#define ATH_MSG_VERBOSE(x)
Definition: AthMsgStreamMacros.h:28
DerivationFramework::JetTrackParticleThinning::finalize
virtual StatusCode finalize() override
Definition: JetTrackParticleThinning.cxx:53
python.utils.AtlRunQueryLookup.mask
string mask
Definition: AtlRunQueryLookup.py:460
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:68
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
lumiFormat.i
int i
Definition: lumiFormat.py:92
beamspotman.n
n
Definition: beamspotman.py:731
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
DerivationFramework::JetTrackParticleThinning::doThinning
virtual StatusCode doThinning() const override
Definition: JetTrackParticleThinning.cxx:61
DerivationFramework::JetTrackParticleThinning::initialize
virtual StatusCode initialize() override
Definition: JetTrackParticleThinning.cxx:32
SG::ReadHandle::isValid
virtual bool isValid() override final
Can the handle be successfully dereferenced?
DerivationFramework::kTrackThinning
@ kTrackThinning
Definition: JetTrackParticleThinning.h:27
JetTrackParticleThinning.h
DerivationFramework::kJetSelection
@ kJetSelection
Definition: JetTrackParticleThinning.h:27
DeMoScan.index
string index
Definition: DeMoScan.py:362
JetContainer.h
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
entries
double entries
Definition: listroot.cxx:49
DerivationFramework::JetTrackParticleThinning::JetTrackParticleThinning
JetTrackParticleThinning(const std::string &t, const std::string &n, const IInterface *p)
Definition: JetTrackParticleThinning.cxx:21
xAOD::JetAttribute::GhostTrack
@ GhostTrack
Definition: JetAttributes.h:252
dqBeamSpot.nEntries
int nEntries
Definition: dqBeamSpot.py:73
DataVector::size
size_type size() const noexcept
Returns the number of elements in the collection.
TrackParticleContainer.h
DerivationFramework::JetTrackParticleThinning::~JetTrackParticleThinning
virtual ~JetTrackParticleThinning()