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

Class for selecting jets that pass cleaning cuts. More...

#include <JetCleaningTool.h>

Inheritance diagram for JetCleaningTool:
Collaboration diagram for JetCleaningTool:

Public Types

enum  CleaningLevel {
  SuperLooseBadLLP , VeryLooseBadLLP , LooseBad , LooseBadLLP ,
  LooseBadTrigger , TightBad , UnknownCut
}
 Levels of cut. More...

Public Member Functions

 JetCleaningTool (const std::string &name="JetCleaningTool")
 Standard constructor.
 JetCleaningTool (const CleaningLevel alevel, const bool doUgly=false)
 Cut-based constructor.
 JetCleaningTool (const std::string &name, const CleaningLevel alevel, const bool doUgly=false)
 Cut and string based constructor.
virtual ~JetCleaningTool ()
 Standard destructor.
virtual StatusCode initialize () override
 Initialize method.
const asg::AcceptInfogetAcceptInfo () const
asg::AcceptData accept (const int isJetClean, const int fmaxIndex) const
 The DFCommonJets decoration accept method.
asg::AcceptData accept (const int isJetClean, const double sumpttrk, const double fmax, const double eta, const double pt, const int fmaxIndex) const
 The DFCommonJets decoration + tight method.
asg::AcceptData accept (const double emf, const double hecf, const double larq, const double hecq, const double sumpttrk, const double eta, const double pt, const double fmax, const double negE, const double AverageLArQF, const int fMaxIndex) const
 The main accept method: the actual cuts are applied here.
asg::AcceptData accept (const xAOD::Jet &jet) const
 The D3PDReader accept method.
int keep (const xAOD::Jet &jet) const final
 Method to select.
virtual StatusCode decorate (const xAOD::JetContainer &jets) const override
 Decorate a jet collection without otherwise modifying it.
bool containsHotCells (const xAOD::Jet &jet, const unsigned int runNumber) const
 Hot cell checks.
CleaningLevel getCutLevel (const std::string &) const
 Helpers for cut names.
std::string getCutName (const CleaningLevel) const
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 StatusCode modify (xAOD::JetContainer &jets) const override final
 Concrete implementation of the function inherited from IJetModifier.
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

StatusCode readHotCells ()
 Hot cells reading helper.
void missingVariable (const std::string &varName) const
Gaudi::Details::PropertyBase & declareGaudiProperty (Gaudi::Property< T, V, H > &hndl, const SG::VarHandleKeyType &)
 specialization for handling Gaudi::Property<SG::VarHandleKey>

Private Attributes

Gaudi::Property< std::string > m_cutName {this, "CutLevel", "" }
 Name of the cut.
CleaningLevel m_cutLevel {LooseBad}
Gaudi::Property< bool > m_doUgly {this, "DoUgly", false}
Gaudi::Property< bool > m_useDecorations {this, "UseDecorations", true}
Gaudi::Property< bool > m_useLooseDecorForTightCut {this, "UseLooseDecorForTightCut", false}
SG::ConstAccessor< char > m_acc_jetClean {"DFCommonJets_jetClean_LooseBad"}
Gaudi::Property< std::string > m_jetContainerName {this, "JetContainer", "", "SG key for input jet container"}
SG::WriteDecorHandleKey< xAOD::JetContainerm_jetCleanKey {this, "JetCleaningName", "isClean", "SG key for output jet cleaning decoration"}
asg::AcceptInfo m_accept
Gaudi::Property< std::string > m_hotCellsFile {this, "HotCellsFile", ""}
 Hot cells caching.
std::unordered_map< unsigned int, std::vector< std::unique_ptr< JCT::HotCell > > > m_hotCellsMap
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

Class for selecting jets that pass cleaning cuts.

Author
Zach Marshall
Date
Feb 2014

Definition at line 39 of file JetCleaningTool.h.

Member Typedef Documentation

◆ StoreGateSvc_t

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

Definition at line 388 of file AthCommonDataStore.h.

Member Enumeration Documentation

◆ CleaningLevel

Levels of cut.

Enumerator
SuperLooseBadLLP 
VeryLooseBadLLP 
LooseBad 
LooseBadLLP 
LooseBadTrigger 
TightBad 
UnknownCut 

Definition at line 46 of file JetCleaningTool.h.

Constructor & Destructor Documentation

◆ JetCleaningTool() [1/3]

JetCleaningTool::JetCleaningTool ( const std::string & name = "JetCleaningTool")

Standard constructor.

Definition at line 89 of file JetCleaningTool.cxx.

90 : asg::AsgTool(name) {}

◆ JetCleaningTool() [2/3]

JetCleaningTool::JetCleaningTool ( const CleaningLevel alevel,
const bool doUgly = false )

Cut-based constructor.

Definition at line 92 of file JetCleaningTool.cxx.

93 : JetCleaningTool( "JetCleaningTool_"+getCutName(alevel) )
94{
95 m_cutLevel=alevel;
96 m_doUgly = doUgly;
97}
JetCleaningTool(const std::string &name="JetCleaningTool")
Standard constructor.
std::string getCutName(const CleaningLevel) const
CleaningLevel m_cutLevel
Gaudi::Property< bool > m_doUgly

◆ JetCleaningTool() [3/3]

JetCleaningTool::JetCleaningTool ( const std::string & name,
const CleaningLevel alevel,
const bool doUgly = false )

Cut and string based constructor.

Definition at line 100 of file JetCleaningTool.cxx.

101 : JetCleaningTool(name)
102{
103 m_cutLevel=alevel;
104 m_doUgly = doUgly;
105}

◆ ~JetCleaningTool()

JetCleaningTool::~JetCleaningTool ( )
virtualdefault

Standard destructor.

Member Function Documentation

◆ accept() [1/4]

asg::AcceptData JetCleaningTool::accept ( const double emf,
const double hecf,
const double larq,
const double hecq,
const double sumpttrk,
const double eta,
const double pt,
const double fmax,
const double negE,
const double AverageLArQF,
const int fMaxIndex ) const

The main accept method: the actual cuts are applied here.

Definition at line 209 of file JetCleaningTool.cxx.

222{
223 asg::AcceptData acceptData (&m_accept);
224 acceptData.setCutResult( "Cleaning", false );
225
226 // -----------------------------------------------------------
227 // Do the actual selection
228 if (pt<DBL_MIN) return acceptData;
229 const double chf=sumpttrk/pt;
230
231 //=============================================================
232 //Run-II ugly cuts
233 //=============================================================
234 if(m_doUgly && fmaxIndex==17) return acceptData;
235
236 //=============================================================
237 //Run-II very loose LLP cuts
238 // From https://indico.cern.ch/event/642438/contributions/2704590/attachments/1514445/2362870/082117a_HCW_NCB_LLP.pdf
239 //=============================================================
241 if (fmax>0.80) return acceptData;
242 if (emf>0.96) return acceptData;
243 acceptData.setCutResult( "Cleaning", true );
244 return acceptData;
245 }
246
247 //=============================================================
248 //Run-II loose cuts
249 //=============================================================
250 //Non-collision background & cosmics
251 const bool useLLP = (LooseBadLLP == m_cutLevel); // LLP cleaning cannot use emf...
252 const bool isTrigger = (LooseBadTrigger == m_cutLevel); // trigger cannot use chf
253 const bool useSuperLLP = (SuperLooseBadLLP == m_cutLevel); // other LLP cleaning...
254 if (!useLLP && !useSuperLLP) {
255 if(!isTrigger && emf<0.05 && chf<0.05 && std::fabs(eta)<2) return acceptData;
256 if(emf<0.05 && std::fabs(eta)>=2) return acceptData;
257 }
258 if(fmax>0.99 && std::fabs(eta)<2) return acceptData;
259 //HEC spike-- gone as of 2017!
260 if(hecf>0.5 && std::fabs(hecq)>0.5 && AverageLArQF/65535>0.8) return acceptData;
261 //EM calo noise
262 if(emf>0.95 && std::fabs(larq)>0.8 && std::fabs(eta)<2.8 && AverageLArQF/65535>0.8) return acceptData;
263 // LLP cleaning uses negative energy cut
264 // (https://indico.cern.ch/event/472320/contribution/8/attachments/1220731/1784456/JetTriggerMeeting_20160102.pdf)
265 if (useLLP && std::fabs(negE*0.001)>4 && fmax >0.85) return acceptData;
266 // another LLP cleaning cutting softer on negative energy
267 if (useSuperLLP && std::fabs(negE*0.001)>60) return acceptData;
268
270 acceptData.setCutResult( "Cleaning", true );
271 return acceptData;
272 }
273
274 //=============================================================
275 //Run-II tight cuts
276 //=============================================================
277 // NCB monojet-style cut in central, EMF cut in forward
278 if (fmax<DBL_MIN) return acceptData;
279 if(std::fabs(eta)<2.4 && chf/fmax<0.1) return acceptData;
280 //if(std::fabs(eta)>=2.4 && emf<0.1) return acceptData;
281 if(TightBad==m_cutLevel){
282 acceptData.setCutResult( "Cleaning", true );
283 return acceptData;
284 }
285
286 // We should never arrive here!
287 ATH_MSG_ERROR( "Unknown cut name: " << getCutName( m_cutLevel ) << " in JetCleaningTool" );
288 return acceptData;
289}
Scalar eta() const
pseudorapidity method
#define ATH_MSG_ERROR(x)
asg::AcceptInfo m_accept

◆ accept() [2/4]

asg::AcceptData JetCleaningTool::accept ( const int isJetClean,
const double sumpttrk,
const double fmax,
const double eta,
const double pt,
const int fmaxIndex ) const

The DFCommonJets decoration + tight method.

Definition at line 172 of file JetCleaningTool.cxx.

179{
180 asg::AcceptData acceptData (&m_accept);
181 acceptData.setCutResult( "Cleaning", false );
182 const double chf=sumpttrk/pt;
183
184 //=============================================================
185 //Run-II ugly cuts
186 //=============================================================
187 if(m_doUgly && fmaxIndex==17) return acceptData;
188
189 //=============================================================
190 //Tight cleaning taken from decoration
191 //=============================================================
192 if(isJetClean==0) return acceptData; //fails Loose cleaning
193 else if (fmax<DBL_MIN) return acceptData;
194 else if(std::fabs(eta)<2.4 && chf/fmax<0.1) return acceptData;
195 else{
196 acceptData.setCutResult( "Cleaning", true );
197 return acceptData;
198 }
199
200 // We should never arrive here!
201 ATH_MSG_ERROR( "Unknown cut name: " << getCutName( m_cutLevel ) << " in JetCleaningTool" );
202 return acceptData;
203
204}

◆ accept() [3/4]

asg::AcceptData JetCleaningTool::accept ( const int isJetClean,
const int fmaxIndex ) const

The DFCommonJets decoration accept method.

Definition at line 145 of file JetCleaningTool.cxx.

146{
147 asg::AcceptData acceptData (&m_accept);
148 acceptData.setCutResult( "Cleaning", false );
149
150 //=============================================================
151 //Run-II ugly cuts
152 //=============================================================
153 if(m_doUgly && fmaxIndex==17) return acceptData;
154
155 //=============================================================
156 //Loose/tight cleaning taken from decoration
157 //=============================================================
158 if(isJetClean==0) return acceptData;
159 else{
160 acceptData.setCutResult( "Cleaning", true );
161 return acceptData;
162 }
163
164 // We should never arrive here!
165 ATH_MSG_ERROR( "Unknown cut name: " << getCutName( m_cutLevel ) << " in JetCleaningTool" );
166 return acceptData;
167
168}

◆ accept() [4/4]

asg::AcceptData JetCleaningTool::accept ( const xAOD::Jet & jet) const

The D3PDReader accept method.

Definition at line 298 of file JetCleaningTool.cxx.

299{
300 std::vector<float> sumPtTrkvec;
302 double sumpttrk = 0;
303 if( ! sumPtTrkvec.empty() ) sumpttrk = sumPtTrkvec[0];
304 // fmax index is not necessarily required
305 // This is only used if doUgly is set
306 // Handle it gracefully if the variable is not present but doUgly is false
307 int FracSamplingMaxIndex = -1;
308 if (!jet.getAttribute(xAOD::JetAttribute::FracSamplingMaxIndex,FracSamplingMaxIndex) && m_doUgly)
309 missingVariable("FracSamplingMaxIndex");
310 // get tight cleaning variables
311 float FracSamplingMax = 0;
312 if (!jet.getAttribute(xAOD::JetAttribute::FracSamplingMax,FracSamplingMax))
313 missingVariable("FracSamplingMax");
314
315 //start jet cleaning
316 int isJetClean = 0;
317 if( m_useDecorations && m_acc_jetClean.isAvailable(jet) ) { //decoration is already available for all jets
318 isJetClean = m_acc_jetClean(jet);
319
320 // compute tight cleaning on top of the loose cleaning decoration
322 return accept (isJetClean, sumpttrk, FracSamplingMax, jet.eta(), jet.pt(), FracSamplingMaxIndex);
323 }
324 else {
325 return accept (isJetClean, FracSamplingMaxIndex);
326 }
327 }
328 else{ //running over AOD, need to use all variables
329 ATH_MSG_DEBUG("DFCommon jet cleaning variable not available ... Using jet cleaning tool");
330 // Get all of the required variables
331 // Do it this way so we can gracefully handle missing variables (rather than segfaults)
332
333 float EMFrac = 0;
335 missingVariable("EMFrac");
336
337 float HECFrac = 0;
339 missingVariable("HECFrac");
340
341 float LArQuality = 0;
342 if (!jet.getAttribute(xAOD::JetAttribute::LArQuality,LArQuality))
343 missingVariable("LArQuality");
344
345
346 float HECQuality = 0;
347 if (!jet.getAttribute(xAOD::JetAttribute::HECQuality,HECQuality))
348 missingVariable("HECQuality");
349
350
351 float NegativeE = 0;
353 missingVariable("NegativeE");
354
355
356
357 float AverageLArQF = 0;
358 if (!jet.getAttribute(xAOD::JetAttribute::AverageLArQF,AverageLArQF))
359 missingVariable("AverageLArQF");
360
361 return accept (EMFrac,
362 HECFrac,
363 LArQuality,
364 HECQuality,
365 sumpttrk,
366 jet.eta(),
367 jet.pt(),
368 FracSamplingMax,
369 NegativeE,
370 AverageLArQF,
371 FracSamplingMaxIndex);}
372}
#define ATH_MSG_DEBUG(x)
asg::AcceptData accept(const int isJetClean, const int fmaxIndex) const
The DFCommonJets decoration accept method.
SG::ConstAccessor< char > m_acc_jetClean
void missingVariable(const std::string &varName) const
Gaudi::Property< bool > m_useDecorations
Gaudi::Property< bool > m_useLooseDecorForTightCut
virtual double pt() const
The transverse momentum ( ) of the particle.
Definition Jet_v1.cxx:44
bool getAttribute(AttributeID type, T &value) const
Retrieve attribute moment by enum.
virtual double eta() const
The pseudorapidity ( ) of the particle.
Definition Jet_v1.cxx:49

◆ containsHotCells()

bool JetCleaningTool::containsHotCells ( const xAOD::Jet & jet,
const unsigned int runNumber ) const

Hot cell checks.

Definition at line 389 of file JetCleaningTool.cxx.

390{
391 // Check if the runNumber contains bad cells
392 std::unordered_map<unsigned int, std::vector<std::unique_ptr<JCT::HotCell>>>::const_iterator hotCells = m_hotCellsMap.find(runNumber);
393 if (hotCells != m_hotCellsMap.end())
394 {
395 // The run contains hot cells
396 // Check if the jet is affected by one of the hot cells
397 for (const std::unique_ptr<JCT::HotCell>& cell : hotCells->second)
398 if (cell->jetAffectedByHotCell(jet))
399 return true;
400 }
401 return false;
402}
std::unordered_map< unsigned int, std::vector< std::unique_ptr< JCT::HotCell > > > m_hotCellsMap

◆ 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 }
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)

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

◆ decorate()

StatusCode JetCleaningTool::decorate ( const xAOD::JetContainer & jets) const
overridevirtual

Decorate a jet collection without otherwise modifying it.

Implements IJetDecorator.

Definition at line 374 of file JetCleaningTool.cxx.

375{
376 ATH_MSG_DEBUG(" Decorating jets with jet cleaning decoration : " << m_jetCleanKey.key());
377
378 SG::WriteDecorHandle<xAOD::JetContainer, bool> cleanHandle(m_jetCleanKey);
379
380 for (const xAOD::Jet *jet : jets) {
381 cleanHandle(*jet) = accept(*jet).getCutResult("Cleaning");
382 }
383
384 return StatusCode::SUCCESS;
385
386}
SG::WriteDecorHandleKey< xAOD::JetContainer > m_jetCleanKey
bool getCutResult(const std::string &cutName) const
Get the result of a cut, based on the cut name (safer)
Definition AcceptData.h:98
Jet_v1 Jet
Definition of the current "jet version".

◆ 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

◆ getAcceptInfo()

const asg::AcceptInfo & JetCleaningTool::getAcceptInfo ( ) const
inline

Definition at line 72 of file JetCleaningTool.h.

73 {
74 return m_accept;
75 }

◆ getCutLevel()

JetCleaningTool::CleaningLevel JetCleaningTool::getCutLevel ( const std::string & s) const

Helpers for cut names.

Definition at line 405 of file JetCleaningTool.cxx.

406{
407 if (s=="SuperLooseBadLLP") return SuperLooseBadLLP;
408 if (s=="VeryLooseBadLLP") return VeryLooseBadLLP;
409 if (s=="LooseBad") return LooseBad;
410 if (s=="LooseBadLLP") return LooseBadLLP;
411 if (s=="LooseBadTrigger") return LooseBadTrigger;
412 if (s=="TightBad") return TightBad;
413 ATH_MSG_ERROR( "Unknown cut level requested: " << s );
414 return UnknownCut;
415}

◆ getCutName()

std::string JetCleaningTool::getCutName ( const CleaningLevel c) const

Definition at line 417 of file JetCleaningTool.cxx.

418{
419 if (c==SuperLooseBadLLP) return "SuperLooseBadLLP";
420 if (c==VeryLooseBadLLP) return "VeryLooseBadLLP";
421 if (c==LooseBad) return "LooseBad";
422 if (c==LooseBadLLP) return "LooseBadLLP";
423 if (c==LooseBadTrigger) return "LooseBadTrigger";
424 if (c==TightBad) return "TightBad";
425 return "UnknownCut";
426}

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

Initialize method.

Reimplemented from asg::AsgTool.

Definition at line 116 of file JetCleaningTool.cxx.

117{
119 ATH_MSG_ERROR( "Tool initialized with unknown cleaning level." );
120 return StatusCode::FAILURE;
121 }
122
124 ATH_MSG_INFO( "Configured with cut level " << getCutName( m_cutLevel ) );
125 std::string jetCleanDFName = "DFCommonJets_jetClean_"+getCutName(m_cutLevel);
126 // if UseLooseDecorForTightCut=true, retrieve loose cleaning decoration to compute tight cleaning
127 if (m_useLooseDecorForTightCut) jetCleanDFName = "DFCommonJets_jetClean_"+getCutName(LooseBad);
128 m_acc_jetClean = jetCleanDFName;
130
131 ATH_MSG_DEBUG( "Initialized decorator name: " << jetCleanDFName );
132
133 m_accept.addCut( "Cleaning", "Cleaning of the jet" );
134
135 // Read in the map of runNumbers to bad cells
137
138 ATH_CHECK(m_jetCleanKey.initialize(!m_jetContainerName.empty()));
139
140 return StatusCode::SUCCESS;
141}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_INFO(x)
Gaudi::Property< std::string > m_jetContainerName
CleaningLevel getCutLevel(const std::string &) const
Helpers for cut names.
Gaudi::Property< std::string > m_cutName
Name of the cut.
StatusCode readHotCells()
Hot cells reading helper.

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

◆ keep()

int JetCleaningTool::keep ( const xAOD::Jet & jet) const
inlinefinalvirtual

Method to select.

Returns true if jet is selected.

Implements IJetSelector.

Definition at line 104 of file JetCleaningTool.h.

105 { return static_cast<bool>(accept(jet)); }

◆ missingVariable()

void JetCleaningTool::missingVariable ( const std::string & varName) const
private

Definition at line 292 of file JetCleaningTool.cxx.

293{
294 ATH_MSG_FATAL(Form("JetCleaningTool failed to retrieve a required variable - please confirm that the xAOD::Jet being passed contains the variable named %s",varName.c_str()));
295 throw std::runtime_error(Form("JetCleaningTool failed to retrieve a required variable - please confirm that the xAOD::Jet being passed contains the variable named %s",varName.c_str()));
296}
#define ATH_MSG_FATAL(x)
str varName
end cluster ToT and charge

◆ modify()

virtual StatusCode IJetDecorator::modify ( xAOD::JetContainer & jets) const
inlinefinaloverridevirtualinherited

Concrete implementation of the function inherited from IJetModifier.

Implements IJetModifier.

Definition at line 32 of file IJetDecorator.h.

32{return decorate(jets);};
virtual StatusCode decorate(const xAOD::JetContainer &jets) const =0
Decorate a jet collection without otherwise modifying it.

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

◆ readHotCells()

StatusCode JetCleaningTool::readHotCells ( )
private

Hot cells reading helper.

Definition at line 430 of file JetCleaningTool.cxx.

431{
432 if (m_hotCellsFile.empty()) return StatusCode::SUCCESS;
433
434 // Ensure that the file exists
436 {
437 ATH_MSG_ERROR("Failed to find hot cells file: " << m_hotCellsFile);
438 return StatusCode::FAILURE;
439 }
440
441 // Now parse the file
442 TEnv readCells;
443 if (readCells.ReadFile(m_hotCellsFile.value().c_str(),kEnvGlobal))
444 {
445 ATH_MSG_ERROR("Cannot read hot cells file: " << m_hotCellsFile);
446 return StatusCode::FAILURE;
447 }
448
449 // Get the list of run numbers
450 const TString runNumbersString = readCells.GetValue("RunNumbers","");
451 if (runNumbersString=="")
452 {
453 ATH_MSG_ERROR("No RunNumbers field was specified in the hot cells file: " << m_hotCellsFile);
454 return StatusCode::FAILURE;
455 }
456
457 // Convert into a vector
458 std::vector<unsigned int> runNumbers = JCT::utils::vectorize<unsigned int>(runNumbersString,", ");
459 if (runNumbers.empty())
460 {
461 ATH_MSG_ERROR("RunNumbers field specified, but value is empty or not unsigned ints for hot cells file: " << m_hotCellsFile);
462 return StatusCode::FAILURE;
463 }
464
465 // Loop over the run numbers
466 for (unsigned int run : runNumbers){
467 std::vector<std::unique_ptr<JCT::HotCell>>& cellVec = m_hotCellsMap[run];
468
469 // The number of hot cells should be below 100 for a given run...
470 for (size_t iCell = 0; iCell < 100; ++iCell)
471 {
472 const TString baseName = Form("Run%u.Cell%zu.",run,iCell);
473
474 // Read the cell info
475 const int layer = readCells.GetValue(baseName+"Layer", -1 );
476 const float minEta = readCells.GetValue(baseName+"EtaMin",-10.);
477 const float maxEta = readCells.GetValue(baseName+"EtaMax", 10.);
478 const float minPhi = readCells.GetValue(baseName+"PhiMin",-10.);
479 const float maxPhi = readCells.GetValue(baseName+"PhiMax", 10.);
480
481 // Check that the input makes sense
482 if (layer < 0 && minEta < -5 && maxEta > 5 && minPhi < -5 && maxPhi > 5)
483 continue;
484 if (layer < 0 || minEta < -5 || maxEta > 5 || minPhi < -5 || maxPhi > 5)
485 {
486 ATH_MSG_ERROR("Partially specified cell - please check the file: " << m_hotCellsFile);
487 ATH_MSG_ERROR(Form("Got Layer=%d, EtaMin=%f, EtaMax=%f, PhiMin=%f, PhiMax=%f",layer,minEta,maxEta,minPhi,maxPhi));
488 return StatusCode::FAILURE;
489 }
490 cellVec.emplace_back(std::make_unique<JCT::HotCell>(layer,minEta,maxEta,minPhi,maxPhi));
491 }
492
493 // Ensure we found the expected run
494 if (cellVec.empty())
495 {
496 ATH_MSG_ERROR("Specified that Run# " << run << " contains hot cells, but did not find any corresponding cells in the file: " << m_hotCellsFile);
497 return StatusCode::FAILURE;
498 }
499 }
500
501 // Done
502 return StatusCode::SUCCESS;
503}
Gaudi::Property< std::string > m_hotCellsFile
Hot cells caching.
bool vectorize(const TString &str, const TString &sep, std::vector< T > &result)
constexpr float maxEta
@ layer
Definition HitInfo.h:79

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

SG::ConstAccessor<char> JetCleaningTool::m_acc_jetClean {"DFCommonJets_jetClean_LooseBad"}
private

Definition at line 124 of file JetCleaningTool.h.

124{"DFCommonJets_jetClean_LooseBad"};

◆ m_accept

asg::AcceptInfo JetCleaningTool::m_accept
private

Definition at line 129 of file JetCleaningTool.h.

◆ m_cutLevel

CleaningLevel JetCleaningTool::m_cutLevel {LooseBad}
private

Definition at line 120 of file JetCleaningTool.h.

120{LooseBad};

◆ m_cutName

Gaudi::Property<std::string> JetCleaningTool::m_cutName {this, "CutLevel", "" }
private

Name of the cut.

Definition at line 119 of file JetCleaningTool.h.

119{this, "CutLevel", "" };

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

Gaudi::Property<bool> JetCleaningTool::m_doUgly {this, "DoUgly", false}
private

Definition at line 121 of file JetCleaningTool.h.

121{this, "DoUgly", false};

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

Gaudi::Property<std::string> JetCleaningTool::m_hotCellsFile {this, "HotCellsFile", ""}
private

Hot cells caching.

Definition at line 132 of file JetCleaningTool.h.

132{this, "HotCellsFile", ""};

◆ m_hotCellsMap

std::unordered_map<unsigned int, std::vector<std::unique_ptr<JCT::HotCell> > > JetCleaningTool::m_hotCellsMap
private

Definition at line 133 of file JetCleaningTool.h.

◆ m_jetCleanKey

SG::WriteDecorHandleKey<xAOD::JetContainer> JetCleaningTool::m_jetCleanKey {this, "JetCleaningName", "isClean", "SG key for output jet cleaning decoration"}
private

Definition at line 128 of file JetCleaningTool.h.

128{this, "JetCleaningName", "isClean", "SG key for output jet cleaning decoration"};

◆ m_jetContainerName

Gaudi::Property<std::string> JetCleaningTool::m_jetContainerName {this, "JetContainer", "", "SG key for input jet container"}
private

Definition at line 127 of file JetCleaningTool.h.

127{this, "JetContainer", "", "SG key for input jet container"};

◆ m_useDecorations

Gaudi::Property<bool> JetCleaningTool::m_useDecorations {this, "UseDecorations", true}
private

Definition at line 122 of file JetCleaningTool.h.

122{this, "UseDecorations", true};

◆ m_useLooseDecorForTightCut

Gaudi::Property<bool> JetCleaningTool::m_useLooseDecorForTightCut {this, "UseLooseDecorForTightCut", false}
private

Definition at line 123 of file JetCleaningTool.h.

123{this, "UseLooseDecorForTightCut", false};

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