ATLAS Offline Software
Loading...
Searching...
No Matches
JetConstituentModSequence Class Reference

#include <JetConstituentModSequence.h>

Inheritance diagram for JetConstituentModSequence:

Public Member Functions

 JetConstituentModSequence (const std::string &name)
StatusCode initialize ()
 Dummy implementation of the initialisation function.
int execute () const
 Method to be called for each event.
virtual void print () const
 Print the state of the tool.
ServiceHandle< StoreGateSvc > & evtStore ()
 The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.
const ServiceHandle< StoreGateSvc > & detStore () const
 The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.
virtual StatusCode sysInitialize () override
 Perform system initialization for an algorithm.
virtual StatusCode sysStart () override
 Handle START transition.
virtual std::vector< Gaudi::DataHandle * > inputHandles () const override
 Return this algorithm's input handles.
virtual std::vector< Gaudi::DataHandle * > outputHandles () const override
 Return this algorithm's output handles.
Gaudi::Details::PropertyBase & declareProperty (Gaudi::Property< T, V, H > &t)
void updateVHKA (Gaudi::Details::PropertyBase &)
MsgStream & msg () const
bool msgLvl (const MSG::Level lvl) const
virtual int inputContainerNames (std::vector< std::string > &connames)
 Method to return the list of input containers.
virtual int outputContainerNames (std::vector< std::string > &connames)
 Method to return the list of output containers.
Additional helper functions, not directly mimicking Athena
template<class T>
const T * getProperty (const std::string &name) const
 Get one of the tool's properties.
const std::string & msg_level_name () const __attribute__((deprecated))
 A deprecated function for getting the message level's name.
const std::string & getName (const void *ptr) const
 Get the name of an object that is / should be in the event store.
SG::sgkey_t getKey (const void *ptr) const
 Get the (hashed) key of an object that is in the event store.

Protected Member Functions

template<class T>
StatusCode copyModRecord (const SG::ReadHandleKey< T > &, const SG::WriteHandleKey< T > &) const
 helper function to cast, shallow copy and record a container.
template<class T, class U>
StatusCode copyModRecordFlowLike (const SG::ReadHandleKey< T > &, const SG::ReadHandleKey< T > &, const SG::WriteHandleKey< T > &, const SG::WriteHandleKey< T > &, const SG::WriteHandleKey< T > &) const
void renounceArray (SG::VarHandleKeyArray &handlesArray)
 remove all handles from I/O resolution
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce (T &h)
void extraDeps_update_handler (Gaudi::Details::PropertyBase &ExtraDeps)
 Add StoreName to extra input/output deps as needed.

Protected Attributes

Gaudi::Property< std::string > m_inputContainer {this, "InputContainer", "", "The input container for the sequence"}
Gaudi::Property< std::string > m_outputContainer = {this, "OutputContainer", "", "The output container for the sequence"}
Gaudi::Property< bool > m_byVertex = {this, "DoByVertex", false, "True if we should match to each primary vertex, not just PV0"}
SG::ReadHandleKey< xAOD::VertexContainerm_vertexContainerKey {this, "VertexContainerKey", "PrimaryVertices", "Reconstructed primary vertex container name"}
unsigned short m_inputType
ToolHandleArray< IJetConstituentModifierm_modifiers {this , "Modifiers" , {} , "List of constit modifier tools."}
ToolHandle< GenericMonitoringToolm_monTool {this,"MonTool","","Monitoring tool"}
bool m_saveAsShallow = true
SG::ReadHandleKey< xAOD::CaloClusterContainerm_inClusterKey {this, "InClusterKey", "", "ReadHandleKey for unmodified CaloClusters"}
SG::WriteHandleKey< xAOD::CaloClusterContainerm_outClusterKey {this, "OutClusterKey", "", "WriteHandleKey for modified CaloClusters"}
SG::ReadHandleKey< xAOD::TrackCaloClusterContainerm_inTCCKey {this, "InTCCKey", "", "ReadHandleKey for unmodified TrackCaloClusters"}
SG::WriteHandleKey< xAOD::TrackCaloClusterContainerm_outTCCKey {this, "OutTCCKey", "", "WriteHandleKey for modified TrackCaloClusters"}
SG::ReadHandleKey< xAOD::PFOContainerm_inChargedPFOKey {this, "InChargedPFOKey", "", "ReadHandleKey for modified Charged PFlow Objects"}
SG::WriteHandleKey< xAOD::PFOContainerm_outChargedPFOKey {this, "OutChargedPFOKey", "", "WriteHandleKey for modified Charged PFlow Objects"}
SG::ReadHandleKey< xAOD::PFOContainerm_inNeutralPFOKey {this, "InNeutralPFOKey", "", "ReadHandleKey for modified Neutral PFlow Objects"}
SG::WriteHandleKey< xAOD::PFOContainerm_outNeutralPFOKey {this, "OutNeutralPFOKey", "", "WriteHandleKey for modified Neutral PFlow Objects"}
SG::WriteHandleKey< xAOD::PFOContainerm_outAllPFOKey {this, "OutAllPFOKey", "", "WriteHandleKey for all modified PFlow Objects"}
SG::ReadHandleKey< xAOD::FlowElementContainerm_inChargedFEKey {this, "InChargedFEKey", "", "ReadHandleKey for modified Charged FlowElements"}
SG::WriteHandleKey< xAOD::FlowElementContainerm_outChargedFEKey {this, "OutChargedFEKey", "", "WriteHandleKey for modified Charged FlowElements"}
SG::ReadHandleKey< xAOD::FlowElementContainerm_inNeutralFEKey {this, "InNeutralFEKey", "", "ReadHandleKey for modified Neutral FlowElements"}
SG::WriteHandleKey< xAOD::FlowElementContainerm_outNeutralFEKey {this, "OutNeutralFEKey", "", "WriteHandleKey for modified Neutral FlowElements"}
SG::WriteHandleKey< xAOD::FlowElementContainerm_outAllFEKey {this, "OutAllFEKey", "", "WriteHandleKey for all modified FlowElements"}
SG::ReadDecorHandleKeyArray< xAOD::FlowElementContainerm_inChargedFEDecorKeys {this, "InChargedFEDecorKeys", {}, "Dummy keys that force all neccessary charged FlowElement decorations to be applied before running this algorithm"}
SG::ReadDecorHandleKeyArray< xAOD::FlowElementContainerm_inNeutralFEDecorKeys {this, "InNeutralFEDecorKeys", {}, "Dummy keys that force all neccessary neutral FlowElement decorations to be applied before running this algorithm"}

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

StoreGateSvc_t m_evtStore
 Pointer to StoreGate (event store by default)
StoreGateSvc_t m_detStore
 Pointer to StoreGate (detector store by default)
std::vector< SG::VarHandleKeyArray * > m_vhka
bool m_varHandleArraysDeclared

Detailed Description

Definition at line 51 of file JetConstituentModSequence.h.

Member Typedef Documentation

◆ StoreGateSvc_t

typedef ServiceHandle<StoreGateSvc> AthCommonDataStore< AthCommonMsg< AlgTool > >::StoreGateSvc_t
privateinherited

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ JetConstituentModSequence()

JetConstituentModSequence::JetConstituentModSequence ( const std::string & name)

Definition at line 26 of file JetConstituentModSequence.cxx.

26 :
27 asg::AsgTool(name) {
28
29#ifdef ASG_TOOL_ATHENA
30 declareInterface<IJetConstituentModifier>(this);
31#endif
32 declareProperty("InputType", m_inputType, "The xAOD type name for the input container.");
33 declareProperty("SaveAsShallow", m_saveAsShallow=true, "Save as shallow copy");
34
35}
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

Member Function Documentation

◆ copyModRecord()

template<class T>
StatusCode JetConstituentModSequence::copyModRecord ( const SG::ReadHandleKey< T > & inKey,
const SG::WriteHandleKey< T > & outKey ) const
protected

helper function to cast, shallow copy and record a container.

Definition at line 121 of file JetConstituentModSequence.h.

122 {
123
124 /* Read in a container of (type is template parameter),
125 optionally modify the elements of this container, and store.
126 This puts a (modified) copy of the container into storegate.
127 */
128
129 auto inHandle = makeHandle(inKey);
130 if(!inHandle.isValid()){
131 ATH_MSG_WARNING("Unable to retrieve input container from " << inKey.key());
132 return StatusCode::FAILURE;
133 }
134
135 std::pair< T*, xAOD::ShallowAuxContainer* > newconstit =
136 xAOD::shallowCopyContainer(*inHandle);
137 newconstit.second->setShallowIO(m_saveAsShallow);
138
139 for (auto t : m_modifiers) {ATH_CHECK(t->process(newconstit.first));}
140
141 auto handle = makeHandle(outKey);
142 ATH_CHECK(handle.record(std::unique_ptr<T>(newconstit.first),
143 std::unique_ptr<xAOD::ShallowAuxContainer>(newconstit.second)));
144
145 xAOD::setOriginalObjectLink(*inHandle, *handle);
146
147 return StatusCode::SUCCESS;
148}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_WARNING(x)
ToolHandleArray< IJetConstituentModifier > m_modifiers
const std::string & key() const
Return the StoreGate ID for the referenced object.
auto makeHandle(const SG::View *view, const KEY &key, const EventContext &ctx)
Create a view handle from a handle key.
Definition ViewHelper.h:273
std::pair< std::unique_ptr< T >, std::unique_ptr< ShallowAuxContainer > > shallowCopyContainer(const T &cont, const EventContext &ctx)
Function making a shallow copy of a constant container.
bool setOriginalObjectLink(const IParticle &original, IParticle &copy)
This function should be used by CP tools when they make a deep copy of an object in their correctedCo...

◆ copyModRecordFlowLike()

template<class T, class U>
StatusCode JetConstituentModSequence::copyModRecordFlowLike ( const SG::ReadHandleKey< T > & inNeutralKey,
const SG::ReadHandleKey< T > & inChargedKey,
const SG::WriteHandleKey< T > & outNeutralKey,
const SG::WriteHandleKey< T > & outChargedKey,
const SG::WriteHandleKey< T > & outAllKey ) const
protected

Definition at line 150 of file JetConstituentModSequence.h.

150 {
151
152 // Cannot be handled the same way as other objects (e.g. clusters),
153 // because the data is split between two containers, but we need
154 // information from both to do the modifications.
155 //
156 // The logic is:
157 // 1. Copy the charged container via a shallow copy
158 // 2. Create N copies of neutral container for by-vertex reconstruction OR
159 // Create a single copy of neutral container for regular jet reconstruction
160 // 3. Merge into a combined view container
161 // 4. Modify the combined container
162
163 // 1. Retrieve the input containers
164 SG::ReadHandle<T> inNeutralHandle = makeHandle(inNeutralKey);
165 SG::ReadHandle<T> inChargedHandle = makeHandle(inChargedKey);
166 if(!inNeutralHandle.isValid()){
167 ATH_MSG_WARNING("Unable to retrieve input containers \""
168 << inNeutralKey.key() << "\" and \""
169 << inChargedKey.key() << "\"");
170 return StatusCode::FAILURE;
171 }
172
173 unsigned numNeutralCopies = 1;
174 if (m_byVertex){
175 // Retrieve Primary Vertices
177 if (!handle.isValid()){
178 ATH_MSG_WARNING(" This event has no primary vertex container" );
179 return StatusCode::FAILURE;
180 }
181
182 const xAOD::VertexContainer* vertices = handle.cptr();
183 if(vertices->empty()){
184 ATH_MSG_WARNING(" Failed to retrieve valid primary vertex container" );
185 return StatusCode::FAILURE;
186 }
187 numNeutralCopies = static_cast<unsigned>(vertices->size());
188 }
189
190 // Copy the input containers individually, set I/O option and record
191 // Charged elements
192 SG::WriteHandle<T> outChargedHandle = makeHandle(outChargedKey);
193
194 std::pair<T*, xAOD::ShallowAuxContainer* > chargedCopy = xAOD::shallowCopyContainer(*inChargedHandle);
195 chargedCopy.second->setShallowIO(m_saveAsShallow);
196 xAOD::setOriginalObjectLink(*inChargedHandle, *chargedCopy.first);
197
198 ATH_CHECK(outChargedHandle.record(std::unique_ptr<T>(chargedCopy.first),
199 std::unique_ptr<xAOD::ShallowAuxContainer>(chargedCopy.second)));
200
201 // Neutral elements
202 SG::WriteHandle<T> outNeutralHandle = makeHandle(outNeutralKey);
203
204 // Shallow copy
205 if (m_saveAsShallow){
206
207 std::pair<T*, xAOD::ShallowAuxContainer* > neutralCopy = xAOD::shallowCopyContainer(*inNeutralHandle);
208 chargedCopy.second->setShallowIO(true);
209 xAOD::setOriginalObjectLink(*inNeutralHandle, *neutralCopy.first);
210
211 ATH_CHECK(outNeutralHandle.record(std::unique_ptr<T>(neutralCopy.first),
212 std::unique_ptr<xAOD::ShallowAuxContainer>(neutralCopy.second)));
213
214 }
215 // Deep copy
216 else{
217 // Define the accessor which will add the index
218 const SG::AuxElement::Accessor<unsigned> copyIndex("ConstituentCopyIndex");
219
220 // Create the new container and its auxiliary store.
221 auto neutralCopies = std::make_unique<T>();
222 auto neutralCopiesAux = std::make_unique<xAOD::AuxContainerBase>();
223 neutralCopies->setStore(neutralCopiesAux.get()); //< Connect the two
224
225 // Create N copies and set copy index
226 // Necessary for by-vertex jet reconstruction
227 for (unsigned i = 0; i < numNeutralCopies; i++){
228 for (const U* fe : *inNeutralHandle) {
229 U* copy = new U();
230 neutralCopies->push_back(copy);
231 *copy = *fe;
232 copyIndex(*copy) = i;
233 xAOD::setOriginalObjectLink(*fe, *copy);
234 }
235 }
236
237 ATH_CHECK(outNeutralHandle.record(std::move(neutralCopies),
238 std::move(neutralCopiesAux))
239 );
240 }
241
242
243 // 2. Set up output handle for merged (view) container and record
244 SG::WriteHandle<T> outAllHandle = makeHandle(outAllKey);
245 ATH_CHECK(outAllHandle.record(std::make_unique<T>(SG::VIEW_ELEMENTS)));
246 (*outAllHandle).assign((*outNeutralHandle).begin(), (*outNeutralHandle).end());
247 (*outAllHandle).insert((*outAllHandle).end(),
248 (*outChargedHandle).begin(),
249 (*outChargedHandle).end());
250
251 // 3. Now process modifications on all elements
252 for (auto t : m_modifiers) {ATH_CHECK(t->process(&*outAllHandle));}
253 return StatusCode::SUCCESS;
254}
size_type size() const noexcept
Returns the number of elements in the collection.
bool empty() const noexcept
Returns true if the collection is empty.
Gaudi::Property< bool > m_byVertex
SG::ReadHandleKey< xAOD::VertexContainer > m_vertexContainerKey
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
virtual bool isValid() override final
Can the handle be successfully dereferenced?
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
@ VIEW_ELEMENTS
this data object is a view, it does not own its elmts
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
bool copy
Definition calibdata.py:26
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".

◆ declareGaudiProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareGaudiProperty ( Gaudi::Property< T, V, H > & hndl,
const SG::VarHandleKeyType &  )
inlineprivateinherited

specialization for handling Gaudi::Property<SG::VarHandleKey>

Definition at line 156 of file AthCommonDataStore.h.

158 {
160 hndl.value(),
161 hndl.documentation());
162
163 }

◆ declareProperty()

Gaudi::Details::PropertyBase & AthCommonDataStore< AthCommonMsg< AlgTool > >::declareProperty ( Gaudi::Property< T, V, H > & t)
inlineinherited

Definition at line 145 of file AthCommonDataStore.h.

145 {
146 typedef typename SG::HandleClassifier<T>::type htype;
148 }
Gaudi::Details::PropertyBase & declareGaudiProperty(Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
specialization for handling Gaudi::Property<SG::VarHandleKey>

◆ detStore()

const ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::detStore ( ) const
inlineinherited

The standard StoreGateSvc/DetectorStore Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 95 of file AthCommonDataStore.h.

◆ evtStore()

ServiceHandle< StoreGateSvc > & AthCommonDataStore< AthCommonMsg< AlgTool > >::evtStore ( )
inlineinherited

The standard StoreGateSvc (event store) Returns (kind of) a pointer to the StoreGateSvc.

Definition at line 85 of file AthCommonDataStore.h.

◆ execute()

int JetConstituentModSequence::execute ( ) const
virtual

Method to be called for each event.

Returns 0 for success.

Implements IJetExecuteTool.

Definition at line 189 of file JetConstituentModSequence.cxx.

189 {
190
191#ifndef XAOD_ANALYSIS
192 // Define monitored quantities
193 auto t_exec = Monitored::Timer<std::chrono::milliseconds>( "TIME_constitmod" );
194#endif
195
196 // Create the shallow copy according to the input type
197 switch(m_inputType){
198
202 if(!sc.isSuccess()) return 1;
203 break;
204 }
205
208 if(!sc.isSuccess()) return 1;
209 break;
210 }
213 if(!sc.isSuccess()) return 1;
214 break;
215 }
216
220 if(!sc.isSuccess()){return 1;}
221 break;
222 }
223
224 default: {
225 ATH_MSG_WARNING( "Unsupported input type " << m_inputType );
226 }
227
228 }
229
230 #ifndef XAOD_ANALYSIS
231 auto mon = Monitored::Group(m_monTool, t_exec);
232 #endif
233 return 0;
234}
static Double_t sc
SG::ReadHandleKey< xAOD::FlowElementContainer > m_inChargedFEKey
StatusCode copyModRecord(const SG::ReadHandleKey< T > &, const SG::WriteHandleKey< T > &) const
helper function to cast, shallow copy and record a container.
ToolHandle< GenericMonitoringTool > m_monTool
SG::WriteHandleKey< xAOD::PFOContainer > m_outChargedPFOKey
SG::WriteHandleKey< xAOD::PFOContainer > m_outNeutralPFOKey
SG::ReadHandleKey< xAOD::PFOContainer > m_inChargedPFOKey
SG::WriteHandleKey< xAOD::PFOContainer > m_outAllPFOKey
SG::WriteHandleKey< xAOD::TrackCaloClusterContainer > m_outTCCKey
SG::WriteHandleKey< xAOD::FlowElementContainer > m_outAllFEKey
SG::WriteHandleKey< xAOD::FlowElementContainer > m_outNeutralFEKey
SG::ReadHandleKey< xAOD::TrackCaloClusterContainer > m_inTCCKey
SG::WriteHandleKey< xAOD::CaloClusterContainer > m_outClusterKey
SG::ReadHandleKey< xAOD::PFOContainer > m_inNeutralPFOKey
StatusCode copyModRecordFlowLike(const SG::ReadHandleKey< T > &, const SG::ReadHandleKey< T > &, const SG::WriteHandleKey< T > &, const SG::WriteHandleKey< T > &, const SG::WriteHandleKey< T > &) const
SG::ReadHandleKey< xAOD::CaloClusterContainer > m_inClusterKey
SG::WriteHandleKey< xAOD::FlowElementContainer > m_outChargedFEKey
SG::ReadHandleKey< xAOD::FlowElementContainer > m_inNeutralFEKey
@ ParticleFlow
The object is a particle-flow object.
Definition ObjectType.h:41
@ FlowElement
The object is a track-calo-cluster.
Definition ObjectType.h:52
@ CaloCluster
The object is a calorimeter cluster.
Definition ObjectType.h:39
@ TrackCaloCluster
The object is a track-calo-cluster.
Definition ObjectType.h:51

◆ extraDeps_update_handler()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::extraDeps_update_handler ( Gaudi::Details::PropertyBase & ExtraDeps)
protectedinherited

Add StoreName to extra input/output deps as needed.

use the logic of the VarHandleKey to parse the DataObjID keys supplied via the ExtraInputs and ExtraOuputs Properties to add the StoreName if it's not explicitly given

◆ getKey()

SG::sgkey_t asg::AsgTool::getKey ( const void * ptr) const
inherited

Get the (hashed) key of an object that is in the event store.

This is a bit of a special one. StoreGateSvc and xAOD::TEvent both provide ways for getting the SG::sgkey_t key for an object that is in the store, based on a bare pointer. But they provide different interfaces for doing so.

In order to allow tools to efficiently perform this operation, they can use this helper function.

See also
asg::AsgTool::getName
Parameters
ptrThe bare pointer to the object that the event store should know about
Returns
The hashed key of the object in the store. If not found, an invalid (zero) key.

Definition at line 119 of file AsgTool.cxx.

119 {
120
121#ifdef XAOD_STANDALONE
122 // In case we use @c xAOD::TEvent, we have a direct function call
123 // for this.
124 return evtStore()->event()->getKey( ptr );
125#else
126 const SG::DataProxy* proxy = evtStore()->proxy( ptr );
127 return ( proxy == nullptr ? 0 : proxy->sgkey() );
128#endif // XAOD_STANDALONE
129 }
ServiceHandle< StoreGateSvc > & evtStore()

◆ getName()

const std::string & asg::AsgTool::getName ( const void * ptr) const
inherited

Get the name of an object that is / should be in the event store.

This is a bit of a special one. StoreGateSvc and xAOD::TEvent both provide ways for getting the std::string name for an object that is in the store, based on a bare pointer. But they provide different interfaces for doing so.

In order to allow tools to efficiently perform this operation, they can use this helper function.

See also
asg::AsgTool::getKey
Parameters
ptrThe bare pointer to the object that the event store should know about
Returns
The string name of the object in the store. If not found, an empty string.

Definition at line 106 of file AsgTool.cxx.

106 {
107
108#ifdef XAOD_STANDALONE
109 // In case we use @c xAOD::TEvent, we have a direct function call
110 // for this.
111 return evtStore()->event()->getName( ptr );
112#else
113 const SG::DataProxy* proxy = evtStore()->proxy( ptr );
114 static const std::string dummy = "";
115 return ( proxy == nullptr ? dummy : proxy->name() );
116#endif // XAOD_STANDALONE
117 }

◆ getProperty()

template<class T>
const T * asg::AsgTool::getProperty ( const std::string & name) const
inherited

Get one of the tool's properties.

◆ initialize()

StatusCode JetConstituentModSequence::initialize ( void )
virtual

Dummy implementation of the initialisation function.

It's here to allow the dual-use tools to skip defining an initialisation function. Since many are doing so...

Reimplemented from asg::AsgTool.

Definition at line 37 of file JetConstituentModSequence.cxx.

37 {
38 ATH_MSG_INFO("Initializing tool " << name() << "...");
39 ATH_MSG_DEBUG("initializing version with data handles");
40
41
42 ATH_CHECK( m_modifiers.retrieve() );
44
45 // Shallow copies are not supported for by-vertex jet reconstruction
46 if (m_byVertex){
47 m_saveAsShallow = false;
48 }
49
50#ifndef XAOD_ANALYSIS
51 ATH_CHECK( m_monTool.retrieve( DisableTool{m_monTool.empty()} ) );
52#endif
53
54 // Set and initialise DataHandleKeys only for the correct input type
55 // Die if the input type is unsupported
56 switch(m_inputType) {
58 {
61
62 ATH_CHECK(m_inClusterKey.initialize());
63 ATH_CHECK(m_outClusterKey.initialize());
64 break;
65 }
67 {
68 std::string inputContainerBase = m_inputContainer;
69 std::string outputContainerBase = m_outputContainer;
70
71 // Know what the user means if they give the full input/output container name in this format
72 size_t pos = inputContainerBase.find("ParticleFlowObjects");
73 if(pos != std::string::npos) inputContainerBase.erase(pos);
74
75 pos = outputContainerBase.find("ParticleFlowObjects");
76 if(pos != std::string::npos) outputContainerBase.erase(pos);
77
78 m_inChargedPFOKey = inputContainerBase + "ChargedParticleFlowObjects";
79 m_inNeutralPFOKey = inputContainerBase + "NeutralParticleFlowObjects";
80
81 m_outChargedPFOKey = outputContainerBase + "ChargedParticleFlowObjects";
82 m_outNeutralPFOKey = outputContainerBase + "NeutralParticleFlowObjects";
83 m_outAllPFOKey = outputContainerBase + "ParticleFlowObjects";
84
85 ATH_CHECK(m_inChargedPFOKey.initialize());
86 ATH_CHECK(m_inNeutralPFOKey.initialize());
87 ATH_CHECK(m_outChargedPFOKey.initialize());
88 ATH_CHECK(m_outNeutralPFOKey.initialize());
89 ATH_CHECK(m_outAllPFOKey.initialize());
90 break;
91 }
92 break;
96
97 ATH_CHECK(m_inTCCKey.initialize());
98 ATH_CHECK(m_outTCCKey.initialize());
99 break;
101 {
102 // TODO: This assumes a PFlow-style neutral and charged collection.
103 // More general FlowElements (e.g. CaloClusters) may necessitate a rework here later.
104
105 const std::string subString = "ParticleFlowObjects";
106 const std::string subStringCharged = "ChargedParticleFlowObjects";
107 const std::string subStringNeutral = "NeutralParticleFlowObjects";
108
109 std::string inputContainerBase = m_inputContainer;
110 std::string outputContainerBase = m_outputContainer;
111
112 m_inChargedFEKey = inputContainerBase;
113 m_inNeutralFEKey = inputContainerBase;
114
115 m_outChargedFEKey = outputContainerBase;
116 m_outNeutralFEKey = outputContainerBase;
117
118 //For both the input and output container basenames we first check if the string
119 //contains "ParticleFlowObjects". If it does we swap this for "ChargedParticleFlowObjects"
120 //and "NeutralParticleFlowObjects" respectively. If it doesn't we just append these two
121 //substrings to the end of the strings.
122
123 size_t pos = inputContainerBase.find(subString);
124 if(pos != std::string::npos) {
125 std::string inChargedString = m_inChargedFEKey.key();
126 m_inChargedFEKey = inChargedString.replace(pos,subString.size(),subStringCharged);
127 std::string inNeutralString = m_inNeutralFEKey.key();
128 m_inNeutralFEKey = inNeutralString.replace(pos,subString.size(),subStringNeutral);
129 }
130 else {
131 m_inChargedFEKey = inputContainerBase + subStringCharged;
132 m_inNeutralFEKey = inputContainerBase + subStringNeutral;
133 }
134
135 pos = outputContainerBase.find(subString);
136 if(pos != std::string::npos) {
137 std::string outChargedString = m_outChargedFEKey.key();
138 m_outChargedFEKey = outChargedString.replace(pos,subString.size(),subStringCharged);
139 std::string outNeutralString = m_outNeutralFEKey.key();
140 m_outNeutralFEKey = outNeutralString.replace(pos,subString.size(),subStringNeutral);
141 }
142 else{
143 m_outChargedFEKey = outputContainerBase + subStringCharged;
144 m_outNeutralFEKey = outputContainerBase + subStringNeutral;
145 }
146
147 //The all FE container is a bit different. If the input container base contains
148 //"ParticleFlowObjects" we do nothing, otherwise we add this string onto the end.
149 pos = outputContainerBase.find(subString);
150 if(pos == std::string::npos) m_outAllFEKey = outputContainerBase + subString;
151 else m_outAllFEKey = outputContainerBase;
152
153 ATH_CHECK(m_inChargedFEKey.initialize());
154 ATH_CHECK(m_inNeutralFEKey.initialize());
155
156 ATH_CHECK(m_outChargedFEKey.initialize());
157 ATH_CHECK(m_outNeutralFEKey.initialize());
158 ATH_CHECK(m_outAllFEKey.initialize());
159
160 // It is enough to initialise the ReadDecorHandleKeys for thread-safety
161 // We are not actually accessing the decorations, only ensuring they are present at runtime
162 // Hence, ReadDecorHandles are not explicitly required to be created
163
164 // Prepend the FE input container name to the decorations
165 for (auto& key : m_inChargedFEDecorKeys) {
166 const std::string keyString = m_inChargedFEKey.key() + "." + key.key();
167 ATH_CHECK(key.assign(keyString));
168
169 }
170
171 for (auto& key : m_inNeutralFEDecorKeys) {
172 const std::string keyString = m_inNeutralFEKey.key() + "." + key.key();
173 ATH_CHECK(key.assign(keyString));
174 }
175
176 ATH_CHECK(m_inChargedFEDecorKeys.initialize());
177 ATH_CHECK(m_inNeutralFEDecorKeys.initialize());
178
179 break;
180 }
181 default:
182 ATH_MSG_ERROR(" Unsupported input type "<< m_inputType );
183 return StatusCode::FAILURE;
184 }
185
186 return StatusCode::SUCCESS;
187}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
Gaudi::Property< std::string > m_outputContainer
SG::ReadDecorHandleKeyArray< xAOD::FlowElementContainer > m_inChargedFEDecorKeys
Gaudi::Property< std::string > m_inputContainer
SG::ReadDecorHandleKeyArray< xAOD::FlowElementContainer > m_inNeutralFEDecorKeys

◆ inputContainerNames()

int IJetExecuteTool::inputContainerNames ( std::vector< std::string > & connames)
virtualinherited

Method to return the list of input containers.

The names of required input containers are appended to connames. Returns nonzero for error. Default returns 0 and adds no names.

Reimplemented in JetRecTool.

Definition at line 11 of file IJetExecuteTool.cxx.

11 {
12 return 0;
13}

◆ inputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::inputHandles ( ) const
overridevirtualinherited

Return this algorithm's input handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ msg()

MsgStream & AthCommonMsg< AlgTool >::msg ( ) const
inlineinherited

Definition at line 24 of file AthCommonMsg.h.

24 {
25 return this->msgStream();
26 }

◆ msg_level_name()

const std::string & asg::AsgTool::msg_level_name ( ) const
inherited

A deprecated function for getting the message level's name.

Instead of using this, weirdly named function, user code should get the string name of the current minimum message level (in case they really need it...), with:

MSG::name( msg().level() )

This function's name doesn't follow the ATLAS coding rules, and as such will be removed in the not too distant future.

Returns
The string name of the current minimum message level that's printed

Definition at line 101 of file AsgTool.cxx.

101 {
102
103 return MSG::name( msg().level() );
104 }
MsgStream & msg() const
const std::string & name(Level lvl)
Convenience function for translating message levels to strings.
Definition MsgLevel.cxx:19

◆ msgLvl()

bool AthCommonMsg< AlgTool >::msgLvl ( const MSG::Level lvl) const
inlineinherited

Definition at line 30 of file AthCommonMsg.h.

30 {
31 return this->msgLevel(lvl);
32 }

◆ outputContainerNames()

int IJetExecuteTool::outputContainerNames ( std::vector< std::string > & connames)
virtualinherited

Method to return the list of output containers.

The names of produced output containers are appended to connames. Returns nonzero for error. Default returns 0 and adds no names.

Reimplemented in JetRecTool.

Definition at line 17 of file IJetExecuteTool.cxx.

17 {
18 return 0;
19}

◆ outputHandles()

virtual std::vector< Gaudi::DataHandle * > AthCommonDataStore< AthCommonMsg< AlgTool > >::outputHandles ( ) const
overridevirtualinherited

Return this algorithm's output handles.

We override this to include handle instances from key arrays if they have not yet been declared. See comments on updateVHKA.

◆ print()

◆ renounce()

std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > AthCommonDataStore< AthCommonMsg< AlgTool > >::renounce ( T & h)
inlineprotectedinherited

Definition at line 380 of file AthCommonDataStore.h.

381 {
382 h.renounce();
384 }
std::enable_if_t< std::is_void_v< std::result_of_t< decltype(&T::renounce)(T)> > &&!std::is_base_of_v< SG::VarHandleKeyArray, T > &&std::is_base_of_v< Gaudi::DataHandle, T >, void > renounce(T &h)

◆ renounceArray()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::renounceArray ( SG::VarHandleKeyArray & handlesArray)
inlineprotectedinherited

remove all handles from I/O resolution

Definition at line 364 of file AthCommonDataStore.h.

364 {
366 }

◆ sysInitialize()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysInitialize ( )
overridevirtualinherited

Perform system initialization for an algorithm.

We override this to declare all the elements of handle key arrays at the end of initialization. See comments on updateVHKA.

Reimplemented in asg::AsgMetadataTool, AthCheckedComponent< AthAlgTool >, AthCheckedComponent<::AthAlgTool >, and DerivationFramework::CfAthAlgTool.

◆ sysStart()

virtual StatusCode AthCommonDataStore< AthCommonMsg< AlgTool > >::sysStart ( )
overridevirtualinherited

Handle START transition.

We override this in order to make sure that conditions handle keys can cache a pointer to the conditions container.

◆ updateVHKA()

void AthCommonDataStore< AthCommonMsg< AlgTool > >::updateVHKA ( Gaudi::Details::PropertyBase & )
inlineinherited

Definition at line 308 of file AthCommonDataStore.h.

308 {
309 // debug() << "updateVHKA for property " << p.name() << " " << p.toString()
310 // << " size: " << m_vhka.size() << endmsg;
311 for (auto &a : m_vhka) {
313 for (auto k : keys) {
314 k->setOwner(this);
315 }
316 }
317 }
std::vector< SG::VarHandleKeyArray * > m_vhka

Member Data Documentation

◆ m_byVertex

Gaudi::Property<bool> JetConstituentModSequence::m_byVertex = {this, "DoByVertex", false, "True if we should match to each primary vertex, not just PV0"}
protected

Definition at line 62 of file JetConstituentModSequence.h.

62{this, "DoByVertex", false, "True if we should match to each primary vertex, not just PV0"};

◆ m_detStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_detStore
privateinherited

Pointer to StoreGate (detector store by default)

Definition at line 393 of file AthCommonDataStore.h.

◆ m_evtStore

StoreGateSvc_t AthCommonDataStore< AthCommonMsg< AlgTool > >::m_evtStore
privateinherited

Pointer to StoreGate (event store by default)

Definition at line 390 of file AthCommonDataStore.h.

◆ m_inChargedFEDecorKeys

SG::ReadDecorHandleKeyArray<xAOD::FlowElementContainer> JetConstituentModSequence::m_inChargedFEDecorKeys {this, "InChargedFEDecorKeys", {}, "Dummy keys that force all neccessary charged FlowElement decorations to be applied before running this algorithm"}
protected

Definition at line 104 of file JetConstituentModSequence.h.

104{this, "InChargedFEDecorKeys", {}, "Dummy keys that force all neccessary charged FlowElement decorations to be applied before running this algorithm"};

◆ m_inChargedFEKey

SG::ReadHandleKey<xAOD::FlowElementContainer> JetConstituentModSequence::m_inChargedFEKey {this, "InChargedFEKey", "", "ReadHandleKey for modified Charged FlowElements"}
protected

Definition at line 94 of file JetConstituentModSequence.h.

94{this, "InChargedFEKey", "", "ReadHandleKey for modified Charged FlowElements"};

◆ m_inChargedPFOKey

SG::ReadHandleKey<xAOD::PFOContainer> JetConstituentModSequence::m_inChargedPFOKey {this, "InChargedPFOKey", "", "ReadHandleKey for modified Charged PFlow Objects"}
protected

Definition at line 86 of file JetConstituentModSequence.h.

86{this, "InChargedPFOKey", "", "ReadHandleKey for modified Charged PFlow Objects"};

◆ m_inClusterKey

SG::ReadHandleKey<xAOD::CaloClusterContainer> JetConstituentModSequence::m_inClusterKey {this, "InClusterKey", "", "ReadHandleKey for unmodified CaloClusters"}
protected

Definition at line 80 of file JetConstituentModSequence.h.

80{this, "InClusterKey", "", "ReadHandleKey for unmodified CaloClusters"};

◆ m_inNeutralFEDecorKeys

SG::ReadDecorHandleKeyArray<xAOD::FlowElementContainer> JetConstituentModSequence::m_inNeutralFEDecorKeys {this, "InNeutralFEDecorKeys", {}, "Dummy keys that force all neccessary neutral FlowElement decorations to be applied before running this algorithm"}
protected

Definition at line 105 of file JetConstituentModSequence.h.

105{this, "InNeutralFEDecorKeys", {}, "Dummy keys that force all neccessary neutral FlowElement decorations to be applied before running this algorithm"};

◆ m_inNeutralFEKey

SG::ReadHandleKey<xAOD::FlowElementContainer> JetConstituentModSequence::m_inNeutralFEKey {this, "InNeutralFEKey", "", "ReadHandleKey for modified Neutral FlowElements"}
protected

Definition at line 97 of file JetConstituentModSequence.h.

97{this, "InNeutralFEKey", "", "ReadHandleKey for modified Neutral FlowElements"};

◆ m_inNeutralPFOKey

SG::ReadHandleKey<xAOD::PFOContainer> JetConstituentModSequence::m_inNeutralPFOKey {this, "InNeutralPFOKey", "", "ReadHandleKey for modified Neutral PFlow Objects"}
protected

Definition at line 89 of file JetConstituentModSequence.h.

89{this, "InNeutralPFOKey", "", "ReadHandleKey for modified Neutral PFlow Objects"};

◆ m_inputContainer

Gaudi::Property<std::string> JetConstituentModSequence::m_inputContainer {this, "InputContainer", "", "The input container for the sequence"}
protected

Definition at line 60 of file JetConstituentModSequence.h.

60{this, "InputContainer", "", "The input container for the sequence"};

◆ m_inputType

unsigned short JetConstituentModSequence::m_inputType
protected

Definition at line 68 of file JetConstituentModSequence.h.

◆ m_inTCCKey

SG::ReadHandleKey<xAOD::TrackCaloClusterContainer> JetConstituentModSequence::m_inTCCKey {this, "InTCCKey", "", "ReadHandleKey for unmodified TrackCaloClusters"}
protected

Definition at line 83 of file JetConstituentModSequence.h.

83{this, "InTCCKey", "", "ReadHandleKey for unmodified TrackCaloClusters"};

◆ m_modifiers

ToolHandleArray<IJetConstituentModifier> JetConstituentModSequence::m_modifiers {this , "Modifiers" , {} , "List of constit modifier tools."}
protected

Definition at line 71 of file JetConstituentModSequence.h.

71{this , "Modifiers" , {} , "List of constit modifier tools."};

◆ m_monTool

ToolHandle<GenericMonitoringTool> JetConstituentModSequence::m_monTool {this,"MonTool","","Monitoring tool"}
protected

Definition at line 74 of file JetConstituentModSequence.h.

74{this,"MonTool","","Monitoring tool"};

◆ m_outAllFEKey

SG::WriteHandleKey<xAOD::FlowElementContainer> JetConstituentModSequence::m_outAllFEKey {this, "OutAllFEKey", "", "WriteHandleKey for all modified FlowElements"}
protected

Definition at line 100 of file JetConstituentModSequence.h.

100{this, "OutAllFEKey", "", "WriteHandleKey for all modified FlowElements"};

◆ m_outAllPFOKey

SG::WriteHandleKey<xAOD::PFOContainer> JetConstituentModSequence::m_outAllPFOKey {this, "OutAllPFOKey", "", "WriteHandleKey for all modified PFlow Objects"}
protected

Definition at line 92 of file JetConstituentModSequence.h.

92{this, "OutAllPFOKey", "", "WriteHandleKey for all modified PFlow Objects"};

◆ m_outChargedFEKey

SG::WriteHandleKey<xAOD::FlowElementContainer> JetConstituentModSequence::m_outChargedFEKey {this, "OutChargedFEKey", "", "WriteHandleKey for modified Charged FlowElements"}
protected

Definition at line 95 of file JetConstituentModSequence.h.

95{this, "OutChargedFEKey", "", "WriteHandleKey for modified Charged FlowElements"};

◆ m_outChargedPFOKey

SG::WriteHandleKey<xAOD::PFOContainer> JetConstituentModSequence::m_outChargedPFOKey {this, "OutChargedPFOKey", "", "WriteHandleKey for modified Charged PFlow Objects"}
protected

Definition at line 87 of file JetConstituentModSequence.h.

87{this, "OutChargedPFOKey", "", "WriteHandleKey for modified Charged PFlow Objects"};

◆ m_outClusterKey

SG::WriteHandleKey<xAOD::CaloClusterContainer> JetConstituentModSequence::m_outClusterKey {this, "OutClusterKey", "", "WriteHandleKey for modified CaloClusters"}
protected

Definition at line 81 of file JetConstituentModSequence.h.

81{this, "OutClusterKey", "", "WriteHandleKey for modified CaloClusters"};

◆ m_outNeutralFEKey

SG::WriteHandleKey<xAOD::FlowElementContainer> JetConstituentModSequence::m_outNeutralFEKey {this, "OutNeutralFEKey", "", "WriteHandleKey for modified Neutral FlowElements"}
protected

Definition at line 98 of file JetConstituentModSequence.h.

98{this, "OutNeutralFEKey", "", "WriteHandleKey for modified Neutral FlowElements"};

◆ m_outNeutralPFOKey

SG::WriteHandleKey<xAOD::PFOContainer> JetConstituentModSequence::m_outNeutralPFOKey {this, "OutNeutralPFOKey", "", "WriteHandleKey for modified Neutral PFlow Objects"}
protected

Definition at line 90 of file JetConstituentModSequence.h.

90{this, "OutNeutralPFOKey", "", "WriteHandleKey for modified Neutral PFlow Objects"};

◆ m_outputContainer

Gaudi::Property<std::string> JetConstituentModSequence::m_outputContainer = {this, "OutputContainer", "", "The output container for the sequence"}
protected

Definition at line 61 of file JetConstituentModSequence.h.

61{this, "OutputContainer", "", "The output container for the sequence"};

◆ m_outTCCKey

SG::WriteHandleKey<xAOD::TrackCaloClusterContainer> JetConstituentModSequence::m_outTCCKey {this, "OutTCCKey", "", "WriteHandleKey for modified TrackCaloClusters"}
protected

Definition at line 84 of file JetConstituentModSequence.h.

84{this, "OutTCCKey", "", "WriteHandleKey for modified TrackCaloClusters"};

◆ m_saveAsShallow

bool JetConstituentModSequence::m_saveAsShallow = true
protected

Definition at line 77 of file JetConstituentModSequence.h.

◆ m_varHandleArraysDeclared

bool AthCommonDataStore< AthCommonMsg< AlgTool > >::m_varHandleArraysDeclared
privateinherited

Definition at line 399 of file AthCommonDataStore.h.

◆ m_vertexContainerKey

SG::ReadHandleKey<xAOD::VertexContainer> JetConstituentModSequence::m_vertexContainerKey {this, "VertexContainerKey", "PrimaryVertices", "Reconstructed primary vertex container name"}
protected

Definition at line 63 of file JetConstituentModSequence.h.

63{this, "VertexContainerKey", "PrimaryVertices", "Reconstructed primary vertex container name"};

◆ m_vhka

std::vector<SG::VarHandleKeyArray*> AthCommonDataStore< AthCommonMsg< AlgTool > >::m_vhka
privateinherited

Definition at line 398 of file AthCommonDataStore.h.


The documentation for this class was generated from the following files: