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

Tool to find jets. More...

#include <JetFinder.h>

Inheritance diagram for JetFinder:
Collaboration diagram for JetFinder:

Public Types

typedef std::vector< std::string > NameList
 Type for ghost labels.

Public Member Functions

 JetFinder (const std::string &name)
virtual StatusCode initialize () override
 Dummy implementation of the initialisation function.
virtual int find (const PseudoJetContainer &cont, xAOD::JetContainer &finalJets, xAOD::JetInput::Type contype) const override
 Method to find jets from a vector of pseudojet inputs.
virtual int findNoSave (const PseudoJetContainer &cont, xAOD::JetContainer &finalJets, xAOD::JetInput::Type contype, fastjet::ClusterSequence *&cs) const override
void save (fastjet::ClusterSequence *pcs) const
bool isVariableR () const
virtual void print () const override
 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
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

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.

Private Types

typedef ServiceHandle< StoreGateSvcStoreGateSvc_t

Private Member Functions

int _find (const PseudoJetContainer &cont, xAOD::JetContainer &finalJets, xAOD::JetInput::Type contype, bool doSave, fastjet::ClusterSequence *&) const
Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

SG::ReadHandleKey< xAOD::EventInfom_eventinfokey {"EventInfo"}
SG::ReadHandleKey< fastjet::ClusterSequence > m_cnameRKey
SG::WriteHandleKey< fastjet::ClusterSequence > m_cnameWKey
std::string m_jetalg
float m_jetrad
float m_minrad
float m_massscale
float m_ptmin
float m_ghostarea
int m_ranopt
ToolHandle< IJetFromPseudojetm_bld
fastjet::JetAlgorithm m_fjalg
bool m_isVariableR
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

Tool to find jets.

Parameters: JetAlgorithm - Algorithm: Kt, AntiKt, CamKt JetRadius - Size parameter (or maximum value for variable-R finding) VariableRMinRadius - Minimum radius for variable-R jet finding VariableRMassScale - Mass scale [MeV] for variable-R jet finding PtMin - PT [MeV] threshold for jet finding. GhostArea - Approximate (starting) area (dy x dphi) for ghost finding RandomOption - Option for area random seets (0=fastjet default, 1=run,event) JetBuilder - Tool used to build jets, interface IJetFromPseudojet (code now moved into this class). Jet active area is evaluated if GhostArea > 0. Variable-R jet finding is performed if VariableRMinRadius >= 0 and VariableRMassScale >= 0. Units are MeV assuming these are the input units

Definition at line 46 of file JetFinder.h.

Member Typedef Documentation

◆ NameList

typedef std::vector<std::string> IJetFinder::NameList
inherited

Type for ghost labels.

Definition at line 36 of file IJetFinder.h.

◆ StoreGateSvc_t

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

Definition at line 388 of file AthCommonDataStore.h.

Constructor & Destructor Documentation

◆ JetFinder()

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

Definition at line 37 of file JetFinder.cxx.

38: AsgTool(name),
39 m_bld("JetFromPseudojet",this),
40 m_fjalg(fastjet::undefined_jet_algorithm),
41 m_isVariableR(false) {
42 declareProperty("JetAlgorithm", m_jetalg="AntiKt");
43 declareProperty("JetRadius", m_jetrad =0.4);
44 declareProperty("VariableRMinRadius", m_minrad =-1.0);
45 declareProperty("VariableRMassScale", m_massscale =-1.0);
46 declareProperty("PtMin", m_ptmin =0.0);
47 declareProperty("GhostArea", m_ghostarea =0.0);
48 declareProperty("RandomOption", m_ranopt =1);
49 declareProperty("JetBuilder", m_bld);
50}
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
std::string m_jetalg
Definition JetFinder.h:84
float m_massscale
Definition JetFinder.h:87
int m_ranopt
Definition JetFinder.h:90
float m_ptmin
Definition JetFinder.h:88
float m_ghostarea
Definition JetFinder.h:89
float m_jetrad
Definition JetFinder.h:85
bool m_isVariableR
Definition JetFinder.h:96
float m_minrad
Definition JetFinder.h:86
ToolHandle< IJetFromPseudojet > m_bld
Definition JetFinder.h:92
fastjet::JetAlgorithm m_fjalg
Definition JetFinder.h:95
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58

Member Function Documentation

◆ _find()

int JetFinder::_find ( const PseudoJetContainer & cont,
xAOD::JetContainer & finalJets,
xAOD::JetInput::Type contype,
bool doSave,
fastjet::ClusterSequence *& pcs ) const
private

Definition at line 119 of file JetFinder.cxx.

123 {
124 if ( m_fjalg == fastjet::undefined_jet_algorithm ) {
125 ATH_MSG_ERROR("Jet finder is not properly initiialized.");
126 return 1;
127 }
128
129 const PseudoJetVector* inps = pjContainer.casVectorPseudoJet();
130
131 ATH_MSG_DEBUG("Jet finding input count: " << inps->size());
132 fastjet::JetDefinition jetdef(m_fjalg, m_jetrad);
133#ifndef NO_JET_VARIABLER
134 const VariableRPlugin* pvrp = nullptr;
135 if ( isVariableR() ) {
136 VariableRPlugin::ClusterType vct = VariableRPlugin::AKTLIKE;
137 switch(m_fjalg) {
138 case fastjet::kt_algorithm: vct = VariableRPlugin::KTLIKE; break;
139 case fastjet::antikt_algorithm: vct = VariableRPlugin::AKTLIKE; break;
140 case fastjet::cambridge_algorithm: vct = VariableRPlugin::CALIKE; break;
141 default:
142 ATH_MSG_ERROR("Invalid algorithm type for variable-R finding.");
143 }
144 pvrp = new VariableRPlugin(m_massscale, m_minrad, m_jetrad, vct, false);
145 jetdef = fastjet::JetDefinition(pvrp);
146 jetdef.delete_plugin_when_unused ();
147 }
148#else
149 if ( isVariableR() ) {
150 ATH_MSG_ERROR("Variable-R jet findng is not supported in theis build.");
151 }
152#endif
153 if ( m_ghostarea <= 0 ) {
154 ATH_MSG_DEBUG("Creating input cluster sequence");
155 pcs = new fastjet::ClusterSequence(*inps, jetdef);
156 } else {
157 fastjet::GhostedAreaSpec gspec(5.0, 1, m_ghostarea);
158 std::vector<int> inseeds;
159 if ( m_ranopt == 1 ) {
160 // Use run/event number as random number seeds.
161
162 auto handle = SG::makeHandle(m_eventinfokey);
163 if (!handle.isValid()){
164 ATH_MSG_ERROR("Unable to retrieve event info");
165 return 1;
166 }
167 const xAOD::EventInfo* pevinfo = handle.cptr();
168
169
170 if ( pevinfo != nullptr ) {
171 auto ievt = pevinfo->eventNumber();
172 auto irun = pevinfo->runNumber();
174 // For MC, use the channel and MC event number
175 ievt = pevinfo->mcEventNumber();
176 irun = pevinfo->mcChannelNumber();
177 }
178 inseeds.push_back(ievt);
179 inseeds.push_back(irun);
180 } else {
181 ATH_MSG_ERROR("Unable to retrieve event info");
182 }
183 } // if (m_ranopt==1)
184 ATH_MSG_DEBUG("Active area specs:");
185 ATH_MSG_DEBUG(" Requested ghost area: " << m_ghostarea);
186 ATH_MSG_DEBUG(" Actual ghost area: " << gspec.actual_ghost_area());
187 ATH_MSG_DEBUG(" Max eta: " << gspec.ghost_etamax());
188 ATH_MSG_DEBUG(" # ghosts: " << gspec.n_ghosts());
189 ATH_MSG_DEBUG(" # rapidity bins: " << gspec.nrap());
190 ATH_MSG_DEBUG(" # phi bins: " << gspec.nphi());
191
192 if ( inseeds.size() == 2 ) {
193 ATH_MSG_DEBUG(" Random seeds: " << inseeds[0] << ", "
194 << inseeds[1]);
195 } else {
196 ATH_MSG_WARNING("Random generator size is not 2: " << inseeds.size());
197 ATH_MSG_DEBUG(" Random seeds: ");
198 for ( auto seed : inseeds ) ATH_MSG_DEBUG(" " << seed);
199 }
200 fastjet::AreaDefinition adef(fastjet::active_area_explicit_ghosts, gspec);
201 //fastjet::AreaDefinition adef(fastjet::active_area, gspec);
202 ATH_MSG_DEBUG("Creating input area cluster sequence");
203 pcs = new fastjet::ClusterSequenceArea(*inps, jetdef,
204 adef.with_fixed_seed(inseeds));
205 }
206
207 ATH_MSG_DEBUG("Calling fastjet");
208 PseudoJetVector outs = sorted_by_pt(pcs->inclusive_jets(m_ptmin));
209 ATH_MSG_DEBUG("Found jet count: " << outs.size());
210 // for ( PseudoJetVector::const_iterator ijet=outs.begin(); ijet!=outs.end(); ++ijet ) {
211 for (const auto & pj: outs ) {
212 xAOD::Jet* pjet = m_bld->add(pj, pjContainer, jets, inputtype);
213
214 // transfer the contituents of pseudojet (which are pseudojets)
215 // to constiuents of jet (which are IParticles)
216 // pjContainer.extractConstituents(*pjet, pj);
217
219 pjet->setAlgorithmType(ialg);
221 if ( isVariableR() ) {
222 pjet->setAttribute("VariableRMinRadius", m_minrad);
223 pjet->setAttribute("VariableRMassScale", m_massscale);
224 }
225 pjet->setAttribute("JetGhostArea", m_ghostarea);
226 }
227 ATH_MSG_DEBUG("Reconstructed jet count: " << jets.size() << " clusterseq="<<pcs);
228
229 for(const xAOD::Jet* j : jets){
230 ATH_MSG_DEBUG("jets reconstructed no of constituents "
231 << j->numConstituents());
232 }
233 // offline stores pcs storegate, trigger: passes it back to caller
234 if(doSave){save(pcs);}
235
236 return 0;
237}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
std::vector< fastjet::PseudoJet > PseudoJetVector
bool isVariableR() const
void save(fastjet::ClusterSequence *pcs) const
SG::ReadHandleKey< xAOD::EventInfo > m_eventinfokey
Definition JetFinder.h:79
uint64_t mcEventNumber() const
The MC generator's event number.
bool eventType(EventType type) const
Check for one particular bitmask value.
@ IS_SIMULATION
true: simulation, false: data
uint32_t runNumber() const
The current event's run number.
uint32_t mcChannelNumber() const
The MC generator's channel number.
uint64_t eventNumber() const
The current event's event number.
void setAlgorithmType(JetAlgorithmType::ID a)
Definition Jet_v1.cxx:258
void setAttribute(const std::string &name, const T &v)
void setSizeParameter(float p)
Definition Jet_v1.cxx:257
SG::ReadCondHandle< T > makeHandle(const SG::ReadCondHandleKey< T > &key, const EventContext &ctx=Gaudi::Hive::currentContext())
ID algId(const std::string &n)
Converts a string into a JetAlgorithmType::ID.
ID
//////////////////////////////////////// JetAlgorithmType::ID defines most common physics jet finding...
Jet_v1 Jet
Definition of the current "jet version".
EventInfo_v1 EventInfo
Definition of the latest event info 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.

◆ 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

◆ find()

int JetFinder::find ( const PseudoJetContainer & cont,
xAOD::JetContainer & finalJets,
xAOD::JetInput::Type contype ) const
overridevirtual

Method to find jets from a vector of pseudojet inputs.

The last arguments are the input type for the found jets and the list of ghost constituent labels. Returns 0 for success.

Implements IJetFinder.

Definition at line 98 of file JetFinder.cxx.

100 {
101
102 constexpr bool doSave = true;
103 fastjet::ClusterSequence* pcs{nullptr};
104 return _find(pjContainer, finalJets, inputtype, doSave, pcs);
105}
int _find(const PseudoJetContainer &cont, xAOD::JetContainer &finalJets, xAOD::JetInput::Type contype, bool doSave, fastjet::ClusterSequence *&) const

◆ findNoSave()

int JetFinder::findNoSave ( const PseudoJetContainer & cont,
xAOD::JetContainer & finalJets,
xAOD::JetInput::Type contype,
fastjet::ClusterSequence *& cs ) const
overridevirtual

Implements IJetFinder.

Definition at line 109 of file JetFinder.cxx.

112 {
113
114 constexpr bool doSave = false;
115 return _find(pjContainer, finalJets, inputtype, doSave, pcs);
116}

◆ 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 JetFinder::initialize ( void )
overridevirtual

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 54 of file JetFinder.cxx.

54 {
55 ATH_MSG_DEBUG("Initializing...");
58 if ( m_fjalg == fastjet::undefined_jet_algorithm ) {
59 ATH_MSG_ERROR("Invalid jet algorithm name: " << m_jetalg);
60 ATH_MSG_ERROR("Allowed values are Kt, CamKt, AntiKt, etc.");
61 return StatusCode::FAILURE;
62 }
63 if ( m_jetrad <=0 ) {
64 ATH_MSG_ERROR("Invalid jet size parameter: " << m_jetrad);
65 return StatusCode::FAILURE;
66 }
67
68 fastjet::JetDefinition jetdef(m_fjalg, m_jetrad);
70 fastjet::ClusterSequence cs(empty, jetdef);
71 cs.inclusive_jets(m_ptmin);
72 m_isVariableR = m_minrad>=0.0 && m_massscale>=0.0;
73#ifdef NO_JET_VARIABLER
74 if ( isVariableR() ) {
75 ATH_MSG_ERROR("Variable-R jet findng is not supported in theis build.");
76 }
77#endif
78
79 // Input DataHandles
80 ATH_CHECK( m_eventinfokey.initialize() );
81
82 std::string sdrop = "ToolSvc.";
83 std::string myname = name();
84 std::string::size_type ipos = myname.find(sdrop);
85 if ( ipos != std::string::npos ) myname.replace(ipos, sdrop.size(), "");
86 std::string cname = "ClusterSequence_JetFinder_" + myname;
87
90 ATH_CHECK( m_cnameRKey.initialize() );
91 ATH_CHECK( m_cnameWKey.initialize() );
92
93 return StatusCode::SUCCESS;
94}
#define ATH_CHECK
Evaluate an expression and check for errors.
static const Attributes_t empty
SG::WriteHandleKey< fastjet::ClusterSequence > m_cnameWKey
Definition JetFinder.h:81
SG::ReadHandleKey< fastjet::ClusterSequence > m_cnameRKey
Definition JetFinder.h:80
fastjet::JetAlgorithm fastJetDef(ID id)

◆ 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.

◆ isVariableR()

bool JetFinder::isVariableR ( ) const

Definition at line 264 of file JetFinder.cxx.

264 {
265 return m_isVariableR;
266}

◆ 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 }

◆ 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()

void JetFinder::print ( ) const
overridevirtual

Print the state of the tool.

Reimplemented from asg::AsgTool.

Definition at line 270 of file JetFinder.cxx.

270 {
271 ATH_MSG_INFO(" Jet algorithm: " << m_jetalg);
272 if ( isVariableR() ) {
273 ATH_MSG_INFO(" Variable-R: true");
274 ATH_MSG_INFO(" min radius: " << m_minrad);
275 ATH_MSG_INFO(" max radius: " << m_jetrad);
276 ATH_MSG_INFO(" mass scale: " << m_massscale);
277 } else {
278 ATH_MSG_INFO(" Jet size parameter: " << m_jetrad);
279 }
280 ATH_MSG_INFO(" Jet min pT [MeV]: " << m_ptmin);
281 ATH_MSG_INFO(" Ghost area: " << m_ghostarea);
282 ATH_MSG_INFO(" Output level: " << MSG::name(msg().level()));
283}
#define ATH_MSG_INFO(x)
MsgStream & msg
Definition testRead.cxx:32

◆ 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 }

◆ save()

void JetFinder::save ( fastjet::ClusterSequence * pcs) const

Definition at line 242 of file JetFinder.cxx.

242 {
243
244 auto handle_out = SG::makeHandle(m_cnameWKey);
245 if ( ! handle_out.record( std::unique_ptr<fastjet::ClusterSequence>(pcs)) ) {
246 ATH_MSG_WARNING("Unable to record " << m_cnameWKey.key());
247 } else {
248 ATH_MSG_DEBUG("Recorded " << m_cnameWKey.key() << " " << pcs );
249 }
250 auto handle_in = SG::makeHandle(m_cnameRKey);
251 bool present = false;
252 if ( handle_in.isValid()) {
253 present = true;
254 }
255 ATH_MSG_DEBUG("Check presence: " << present);
256 ATH_MSG_DEBUG("Will delete self: " << pcs->will_delete_self_when_unused());
257 const fastjet::SharedPtr<fastjet::PseudoJetStructureBase>& shrptr = pcs->structure_shared_ptr();
258 ATH_MSG_DEBUG(" Pointer: " << shrptr.get());
259 ATH_MSG_DEBUG(" Use count: " << shrptr.use_count());
260}

◆ 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_bld

ToolHandle<IJetFromPseudojet> JetFinder::m_bld
private

Definition at line 92 of file JetFinder.h.

◆ m_cnameRKey

SG::ReadHandleKey<fastjet::ClusterSequence> JetFinder::m_cnameRKey
private

Definition at line 80 of file JetFinder.h.

◆ m_cnameWKey

SG::WriteHandleKey<fastjet::ClusterSequence> JetFinder::m_cnameWKey
private

Definition at line 81 of file JetFinder.h.

◆ 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_eventinfokey

SG::ReadHandleKey<xAOD::EventInfo> JetFinder::m_eventinfokey {"EventInfo"}
private

Definition at line 79 of file JetFinder.h.

79{"EventInfo"};

◆ 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_fjalg

fastjet::JetAlgorithm JetFinder::m_fjalg
private

Definition at line 95 of file JetFinder.h.

◆ m_ghostarea

float JetFinder::m_ghostarea
private

Definition at line 89 of file JetFinder.h.

◆ m_isVariableR

bool JetFinder::m_isVariableR
private

Definition at line 96 of file JetFinder.h.

◆ m_jetalg

std::string JetFinder::m_jetalg
private

Definition at line 84 of file JetFinder.h.

◆ m_jetrad

float JetFinder::m_jetrad
private

Definition at line 85 of file JetFinder.h.

◆ m_massscale

float JetFinder::m_massscale
private

Definition at line 87 of file JetFinder.h.

◆ m_minrad

float JetFinder::m_minrad
private

Definition at line 86 of file JetFinder.h.

◆ m_ptmin

float JetFinder::m_ptmin
private

Definition at line 88 of file JetFinder.h.

◆ m_ranopt

int JetFinder::m_ranopt
private

Definition at line 90 of file JetFinder.h.

◆ m_varHandleArraysDeclared

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

Definition at line 399 of file AthCommonDataStore.h.

◆ 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: