ATLAS Offline Software
Loading...
Searching...
No Matches
PseudoJetAlgorithm.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5// PseudoJetAlgorithm.cxx
6
12
13// Fixed value by which to scale ghost kinematics
14constexpr float ghostscale = 1e-40;
15
16//**********************************************************************
17
19 ATH_MSG_INFO("Initializing " << name() << " with properties:");
20
21 // This is horrible, but still necessary for the time being
22 // PJG needs to know if this is the basic EMTopo cluster collection
23 // in order to change the cluster signal state.
24 if ( m_label == "EMTopo") m_emtopo = true;
25 // PFlow containers need to have PV matching applied
26 if ( std::string(m_label).find("PFlow") != std::string::npos) m_pflow = true;
27 if ( std::string(m_label).find("UFO") != std::string::npos) m_ufo = true;
28
29 // "Ghost" in output collection name? If so is a ghost collection.
30 m_isGhost = (m_outcoll.key()).find("Ghost") != std::string::npos;
31
32 ATH_MSG_INFO(" " << m_label.value() << ": " << m_incoll.key() << " -> " << m_outcoll.key());
33
34 print();
35
36 if(m_incoll.key().empty() || m_outcoll.key().empty()) {
37 ATH_MSG_ERROR("Either input or output collection is empty!");
38 return StatusCode::FAILURE;
39 }
40
41 ATH_CHECK( m_incoll.initialize() );
42 ATH_CHECK( m_outcoll.initialize() );
43
45
46 return StatusCode::SUCCESS;
47}
48
49
50//**********************************************************************
51
52StatusCode PseudoJetAlgorithm::execute(const EventContext& ctx) const {
53 ATH_MSG_VERBOSE("Executing " << name() << "...");
54
55 // build and record the container
56 auto incoll = SG::makeHandle(m_incoll, ctx);
57 if( !incoll.isValid() ) {
58 // Return SUCCESS to avoid crashing T0 jobs
59 ATH_MSG_WARNING("Failed to retrieve " << m_incoll.key() << " for PseudoJet creation!" );
60 return StatusCode::SUCCESS;
61 }
62 ATH_MSG_DEBUG("Retrieved xAOD container " << m_incoll.key() << " of size " << incoll->size()
63 << ", isGhost=" << m_isGhost);
64
65 ATH_MSG_DEBUG("Creating PseudoJetContainer...");
66 std::unique_ptr<PseudoJetContainer> pjcont( createPJContainer(*incoll, ctx) );
67
68 auto outcoll = SG::makeHandle(m_outcoll, ctx);
69 ATH_MSG_DEBUG("Created new PseudoJetContainer \"" << m_outcoll.key() << "\" with size " << pjcont->size());
70
71 ATH_CHECK( outcoll.record(std::move(pjcont)) );
72
73 return StatusCode::SUCCESS;
74}
75
76
77std::unique_ptr<PseudoJetContainer> PseudoJetAlgorithm::createPJContainer(const xAOD::IParticleContainer& cont, const EventContext& ctx) const {
78 // create PseudoJets
79 std::vector<fastjet::PseudoJet> vpj;
80
81 #ifndef GENERATIONBASE
82 if (m_byVertex){
84 vpj = createPseudoJets(cont, pvs.cptr());
85 }
86 else{
87 vpj = createPseudoJets(cont);
88 }
89 #else
90 [[maybe_unused]] const EventContext& unused_ctx = ctx;
91 vpj = createPseudoJets(cont);
92 #endif
93
94 // create an extractor to attach to the PJContainer -- this will be used by clients
95 auto extractor = std::make_unique<IParticleExtractor>(&cont, m_label, m_isGhost);
96 ATH_MSG_DEBUG("Created extractor: " << extractor->toString(0));
97
98 // ghostify the pseudojets if necessary
99 if(m_isGhost){
100 for(fastjet::PseudoJet& pj : vpj) {pj *= ghostscale;}
101 }
102
103 // Put the PseudoJetContainer together
104 auto pjcont = std::make_unique<PseudoJetContainer>(std::move(extractor), vpj);
105 ATH_MSG_DEBUG("New PseudoJetContainer size " << pjcont->size());
106
107 return pjcont;
108}
109
110
111std::vector<fastjet::PseudoJet>
113#ifndef GENERATIONBASE
116#endif
118}
119
120#ifndef GENERATIONBASE
121std::vector<fastjet::PseudoJet>
125#endif
126
127//**********************************************************************
128
130 std::string sskip = m_skipNegativeEnergy ? "true" : "false";
131 ATH_MSG_DEBUG(" Skip negative E: " << sskip);
132 ATH_MSG_DEBUG(" Is EMTopo: " << m_emtopo);
133 ATH_MSG_DEBUG(" Is PFlow: " << m_pflow);
134 ATH_MSG_DEBUG(" Is UFO: " << m_ufo);
135 ATH_MSG_DEBUG(" Is ghost: " << m_isGhost);
136 ATH_MSG_DEBUG(" Treat negative E as ghost: " << m_negEnergyAsGhosts.value());
137 ATH_MSG_DEBUG(" Running by-vertex reco: " << m_byVertex.value());
138 ATH_MSG_DEBUG(" Vertex input container: " << m_vertexContainer_key.key());
139
140 if(m_pflow){
141 ATH_MSG_DEBUG(" Use charged FEs: " << m_useCharged.value());
142 ATH_MSG_DEBUG(" Use neutral FEs: " << m_useNeutral.value());
143 ATH_MSG_DEBUG("Use charged PV FEs: " << m_useChargedPV.value());
144 ATH_MSG_DEBUG(" PU sideband def: " << m_useChargedPUsideband.value());
145 }
146 if(m_ufo){
147 ATH_MSG_DEBUG(" Use charged UFOs: " << m_useCharged.value());
148 ATH_MSG_DEBUG(" Use neutral UFOs: " << m_useNeutral.value());
149 }
150}
151
152//**********************************************************************
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Handle class for reading from StoreGate.
Handle class for recording to StoreGate.
constexpr float ghostscale
Gaudi::Property< bool > m_useNeutral
Flag to define if neutral PFOs / FEs should be considered.
bool m_pflow
True if inputs are EM-scale topo clusters.
Gaudi::Property< bool > m_useChargedPV
Flag to define if charged PFOs / FEs should be matched to PV.
std::vector< fastjet::PseudoJet > createPseudoJets(const xAOD::IParticleContainer &) const
bool m_ufo
True if inputs are PFlow.
Gaudi::Property< std::string > m_label
Label for the collection.
Gaudi::Property< bool > m_useChargedPUsideband
Flag for PFlow sideband definition.
Gaudi::Property< bool > m_useCharged
Flag to define if charged PFOs / FEs should be considered.
SG::WriteHandleKey< PseudoJetContainer > m_outcoll
Output collection name.
bool m_emtopo
Determines whether the PJs should be made ghosts.
bool m_isGhost
Internal steering flags Set in initialize()
SG::ReadHandleKey< xAOD::VertexContainer > m_vertexContainer_key
Gaudi::Property< bool > m_skipNegativeEnergy
Flag indicating to skip objects with E<0.
virtual std::unique_ptr< PseudoJetContainer > createPJContainer(const xAOD::IParticleContainer &cont, const EventContext &ctx) const
Method to construct the PseudoJetContainer and record in StoreGate.
Gaudi::Property< bool > m_byVertex
Flag for by-vertex jet reconstruction.
virtual void print() const
Dump to properties to the log.
virtual StatusCode execute(const EventContext &ctx) const override final
virtual StatusCode initialize() override final
Athena algorithm's Hooks.
Gaudi::Property< bool > m_negEnergyAsGhosts
Flag indicating to treat objects with E<0 as ghosts (useful for HI)
SG::ReadHandleKey< xAOD::IParticleContainer > m_incoll
Input collection name.
std::string find(const std::string &s)
return a remapped string
Definition hcg.cxx:138
std::vector< fastjet::PseudoJet > EMToposToPJs(const xAOD::IParticleContainer &ips, bool skipNegativeEnergy)
std::vector< fastjet::PseudoJet > PFlowsToPJs(const xAOD::IParticleContainer &ips, bool skipNegativeEnergy, bool useChargedPFOs, bool useNeutralPFOs, bool useChargedPV, bool useChargedPUsideband, bool isUFO)
std::vector< fastjet::PseudoJet > IParticlesToPJs(const xAOD::IParticleContainer &ips, bool skipNegativeEnergy)
std::vector< fastjet::PseudoJet > ByVertexPFlowsToPJs(const xAOD::IParticleContainer &ips, const xAOD::VertexContainer *pvs, bool skipNegativeEnergy, bool useChargedPFOs, bool useNeutralPFOs, bool isUFO)
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
DataVector< IParticle > IParticleContainer
Simple convenience declaration of IParticleContainer.