ATLAS Offline Software
Loading...
Searching...
No Matches
JetCleaningTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5/******************************************************************************
6Name: JetCleaningTool
7
8Author: Zach Marshall
9Created: Feb 2014
10
11Description: Class for selecting jets that pass some cleaning cuts
12******************************************************************************/
13
16
17// This class header and package headers
20
21// xAOD/ASG includes
23
24// STL includes
25#include <iostream>
26#include <cmath>
27#include <cfloat>
28#include <stdexcept>
29#include <utility>
30// ROOT includes
31#include "TEnv.h"
32
33//=============================================================================
34// Constructors
35//=============================================================================
36namespace JCT
37{
39{
40 public:
41 HotCell(const int layer, const float etaMin, const float etaMax, const float phiMin, const float phiMax);
42 virtual ~HotCell() {}
43 bool jetAffectedByHotCell(const xAOD::Jet& jet) const;
44 private:
45 const int m_layer;
46 const float m_etaMin;
47 const float m_etaMax;
48 const float m_phiMin;
49 const float m_phiMax;
51};
52
53HotCell::HotCell(const int layer, const float etaMin, const float etaMax, const float phiMin, const float phiMax)
54 : asg::AsgMessaging("HotCell")
55 , m_layer(layer)
56 , m_etaMin(etaMin)
57 , m_etaMax(etaMax)
58 , m_phiMin(phiMin)
59 , m_phiMax(phiMax) { }
60
62{
63 // First check if the jet points to the cell
64 const float eta = jet.eta();
65 const float phi = jet.phi();
66 if ( (m_etaMin < eta && eta < m_etaMax) && (m_phiMin < phi && phi < m_phiMax) )
67 {
68 // It points to the cell, now check if the maximum layer is the hot cell or if at least 40% of the energy is the hot cell
69 const int fmaxIndex = jet.getAttribute<int>(xAOD::JetAttribute::FracSamplingMaxIndex);
70 float ePerSamp = 0;
71 if( m_ePerSamp.isAvailable(jet) )
72 ePerSamp = m_ePerSamp(jet).at(m_layer)/jet.e();
73 else
74 ATH_MSG_WARNING("Could not retrieve EnergyPerSampling from jet, cleaning performance may be reduced");
75
76 if (fmaxIndex == m_layer || ePerSamp > 0.4)
77 return true;
78 }
79 return false;
80}
81
82} // end JCT namespace
83
84
85
86//=============================================================================
87// Constructors
88//=============================================================================
89JetCleaningTool::JetCleaningTool(const std::string& name)
90 : asg::AsgTool(name) {}
91
92JetCleaningTool::JetCleaningTool(const CleaningLevel alevel, const bool doUgly)
93 : JetCleaningTool( "JetCleaningTool_"+getCutName(alevel) )
94{
95 m_cutLevel=alevel;
96 m_doUgly = doUgly;
97}
98
100JetCleaningTool::JetCleaningTool(const std::string& name , const CleaningLevel alevel, const bool doUgly)
101 : JetCleaningTool(name)
102{
103 m_cutLevel=alevel;
104 m_doUgly = doUgly;
105}
106
107
108//=============================================================================
109// Destructor
110//=============================================================================
112
113//=============================================================================
114// Initialize
115//=============================================================================
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}
142//===============================================================
143// Calculate the accept from the DFCommonJets_jetClean decorator
144//===============================================================
145asg::AcceptData JetCleaningTool::accept( const int isJetClean, const int fmaxIndex ) const
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 acceptData.setCutResult( "Cleaning", true );
160 return acceptData;
161}
162//===============================================================
163// Calculate tight cleaning from loose decoration + variables
164//===============================================================
166 const double sumpttrk, //in MeV, same as sumpttrk
167 const double fmax,
168 const double eta,
169 const double pt,
170 const int fmaxIndex
171 ) const
172{
173 asg::AcceptData acceptData (&m_accept);
174 acceptData.setCutResult( "Cleaning", false );
175 const double chf=sumpttrk/pt;
176
177 //=============================================================
178 //Run-II ugly cuts
179 //=============================================================
180 if(m_doUgly && fmaxIndex==17) return acceptData;
181
182 //=============================================================
183 //Tight cleaning taken from decoration
184 //=============================================================
185 if(isJetClean==0) return acceptData; //fails Loose cleaning
186 else if (fmax<DBL_MIN) return acceptData;
187 else if(std::fabs(eta)<2.4 && chf/fmax<0.1) return acceptData;
188 //
189 acceptData.setCutResult( "Cleaning", true );
190 return acceptData;
191
192}
193
194//=============================================================================
195// Calculate the actual accept of each cut individually.
196//=============================================================================
198 const double hecf,
199 const double larq,
200 const double hecq,
201 //const double time, //in ns
202 const double sumpttrk, //in MeV, same as sumpttrk
203 const double eta, //emscale Eta
204 const double pt, //in MeV, same as sumpttrk
205 const double fmax,
206 const double negE , //in MeV
207 const double AverageLArQF,
208 const int fmaxIndex
209 ) const
210{
211 asg::AcceptData acceptData (&m_accept);
212 acceptData.setCutResult( "Cleaning", false );
213
214 // -----------------------------------------------------------
215 // Do the actual selection
216 if (pt<DBL_MIN) return acceptData;
217 const double chf=sumpttrk/pt;
218
219 //=============================================================
220 //Run-II ugly cuts
221 //=============================================================
222 if(m_doUgly && fmaxIndex==17) return acceptData;
223
224 //=============================================================
225 //Run-II very loose LLP cuts
226 // From https://indico.cern.ch/event/642438/contributions/2704590/attachments/1514445/2362870/082117a_HCW_NCB_LLP.pdf
227 //=============================================================
229 if (fmax>0.80) return acceptData;
230 if (emf>0.96) return acceptData;
231 acceptData.setCutResult( "Cleaning", true );
232 return acceptData;
233 }
234
235 //=============================================================
236 //Run-II loose cuts
237 //=============================================================
238 //Non-collision background & cosmics
239 const bool useLLP = (LooseBadLLP == m_cutLevel); // LLP cleaning cannot use emf...
240 const bool isTrigger = (LooseBadTrigger == m_cutLevel); // trigger cannot use chf
241 const bool useSuperLLP = (SuperLooseBadLLP == m_cutLevel); // other LLP cleaning...
242 if (!useLLP && !useSuperLLP) {
243 if(!isTrigger && emf<0.05 && chf<0.05 && std::fabs(eta)<2) return acceptData;
244 if(emf<0.05 && std::fabs(eta)>=2) return acceptData;
245 }
246 if(fmax>0.99 && std::fabs(eta)<2) return acceptData;
247 //HEC spike-- gone as of 2017!
248 if(hecf>0.5 && std::fabs(hecq)>0.5 && AverageLArQF/65535>0.8) return acceptData;
249 //EM calo noise
250 if(emf>0.95 && std::fabs(larq)>0.8 && std::fabs(eta)<2.8 && AverageLArQF/65535>0.8) return acceptData;
251 // LLP cleaning uses negative energy cut
252 // (https://indico.cern.ch/event/472320/contribution/8/attachments/1220731/1784456/JetTriggerMeeting_20160102.pdf)
253 if (useLLP && std::fabs(negE*0.001)>4 && fmax >0.85) return acceptData;
254 // another LLP cleaning cutting softer on negative energy
255 if (useSuperLLP && std::fabs(negE*0.001)>60) return acceptData;
256
258 acceptData.setCutResult( "Cleaning", true );
259 return acceptData;
260 }
261
262 //=============================================================
263 //Run-II tight cuts
264 //=============================================================
265 // NCB monojet-style cut in central, EMF cut in forward
266 if (fmax<DBL_MIN) return acceptData;
267 if(std::fabs(eta)<2.4 && chf/fmax<0.1) return acceptData;
268 //if(std::fabs(eta)>=2.4 && emf<0.1) return acceptData;
269 if(TightBad==m_cutLevel){
270 acceptData.setCutResult( "Cleaning", true );
271 return acceptData;
272 }
273
274 // We should never arrive here!
275 ATH_MSG_ERROR( "Unknown cut name: " << getCutName( m_cutLevel ) << " in JetCleaningTool" );
276 return acceptData;
277}
278
279
280void JetCleaningTool::missingVariable(const std::string& varName) const
281{
282 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()));
283 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()));
284}
285
287{
288 std::vector<float> sumPtTrkvec;
289 jet.getAttribute( xAOD::JetAttribute::SumPtTrkPt500, sumPtTrkvec );
290 double sumpttrk = 0;
291 if( ! sumPtTrkvec.empty() ) sumpttrk = sumPtTrkvec[0];
292 // fmax index is not necessarily required
293 // This is only used if doUgly is set
294 // Handle it gracefully if the variable is not present but doUgly is false
295 int FracSamplingMaxIndex = -1;
296 if (!jet.getAttribute(xAOD::JetAttribute::FracSamplingMaxIndex,FracSamplingMaxIndex) && m_doUgly)
297 missingVariable("FracSamplingMaxIndex");
298 // get tight cleaning variables
299 float FracSamplingMax = 0;
300 if (!jet.getAttribute(xAOD::JetAttribute::FracSamplingMax,FracSamplingMax))
301 missingVariable("FracSamplingMax");
302
303 //start jet cleaning
304 int isJetClean = 0;
305 if( m_useDecorations && m_acc_jetClean.isAvailable(jet) ) { //decoration is already available for all jets
306 isJetClean = m_acc_jetClean(jet);
307
308 // compute tight cleaning on top of the loose cleaning decoration
310 return accept (isJetClean, sumpttrk, FracSamplingMax, jet.eta(), jet.pt(), FracSamplingMaxIndex);
311 }
312 else {
313 return accept (isJetClean, FracSamplingMaxIndex);
314 }
315 }
316 else{ //running over AOD, need to use all variables
317 ATH_MSG_DEBUG("DFCommon jet cleaning variable not available ... Using jet cleaning tool");
318 // Get all of the required variables
319 // Do it this way so we can gracefully handle missing variables (rather than segfaults)
320
321 float EMFrac = 0;
322 if (!jet.getAttribute(xAOD::JetAttribute::EMFrac,EMFrac))
323 missingVariable("EMFrac");
324
325 float HECFrac = 0;
326 if (!jet.getAttribute(xAOD::JetAttribute::HECFrac,HECFrac))
327 missingVariable("HECFrac");
328
329 float LArQuality = 0;
330 if (!jet.getAttribute(xAOD::JetAttribute::LArQuality,LArQuality))
331 missingVariable("LArQuality");
332
333
334 float HECQuality = 0;
335 if (!jet.getAttribute(xAOD::JetAttribute::HECQuality,HECQuality))
336 missingVariable("HECQuality");
337
338
339 float NegativeE = 0;
340 if (!jet.getAttribute(xAOD::JetAttribute::NegativeE,NegativeE))
341 missingVariable("NegativeE");
342
343
344
345 float AverageLArQF = 0;
346 if (!jet.getAttribute(xAOD::JetAttribute::AverageLArQF,AverageLArQF))
347 missingVariable("AverageLArQF");
348
349 return accept (EMFrac,
350 HECFrac,
351 LArQuality,
352 HECQuality,
353 sumpttrk,
354 jet.eta(),
355 jet.pt(),
356 FracSamplingMax,
357 NegativeE,
358 AverageLArQF,
359 FracSamplingMaxIndex);}
360}
361
363{
364 ATH_MSG_DEBUG(" Decorating jets with jet cleaning decoration : " << m_jetCleanKey.key());
365
367
368 for (const xAOD::Jet *jet : jets) {
369 cleanHandle(*jet) = accept(*jet).getCutResult("Cleaning");
370 }
371
372 return StatusCode::SUCCESS;
373
374}
375
377bool JetCleaningTool::containsHotCells( const xAOD::Jet& jet, const unsigned int runNumber) const
378{
379 // Check if the runNumber contains bad cells
380 std::unordered_map<unsigned int, std::vector<std::unique_ptr<JCT::HotCell>>>::const_iterator hotCells = m_hotCellsMap.find(runNumber);
381 if (hotCells != m_hotCellsMap.end())
382 {
383 // The run contains hot cells
384 // Check if the jet is affected by one of the hot cells
385 for (const std::unique_ptr<JCT::HotCell>& cell : hotCells->second)
386 if (cell->jetAffectedByHotCell(jet))
387 return true;
388 }
389 return false;
390}
391
394{
395 if (s=="SuperLooseBadLLP") return SuperLooseBadLLP;
396 if (s=="VeryLooseBadLLP") return VeryLooseBadLLP;
397 if (s=="LooseBad") return LooseBad;
398 if (s=="LooseBadLLP") return LooseBadLLP;
399 if (s=="LooseBadTrigger") return LooseBadTrigger;
400 if (s=="TightBad") return TightBad;
401 return UnknownCut;
402}
403
405{
406 if (c==SuperLooseBadLLP) return "SuperLooseBadLLP";
407 if (c==VeryLooseBadLLP) return "VeryLooseBadLLP";
408 if (c==LooseBad) return "LooseBad";
409 if (c==LooseBadLLP) return "LooseBadLLP";
410 if (c==LooseBadTrigger) return "LooseBadTrigger";
411 if (c==TightBad) return "TightBad";
412 return "UnknownCut";
413}
414
415
418{
419 if (m_hotCellsFile.empty()) return StatusCode::SUCCESS;
420
421 // Ensure that the file exists
423 {
424 ATH_MSG_ERROR("Failed to find hot cells file: " << m_hotCellsFile);
425 return StatusCode::FAILURE;
426 }
427
428 // Now parse the file
429 TEnv readCells;
430 if (readCells.ReadFile(m_hotCellsFile.value().c_str(),kEnvGlobal))
431 {
432 ATH_MSG_ERROR("Cannot read hot cells file: " << m_hotCellsFile);
433 return StatusCode::FAILURE;
434 }
435
436 // Get the list of run numbers
437 const TString runNumbersString = readCells.GetValue("RunNumbers","");
438 if (runNumbersString=="")
439 {
440 ATH_MSG_ERROR("No RunNumbers field was specified in the hot cells file: " << m_hotCellsFile);
441 return StatusCode::FAILURE;
442 }
443
444 // Convert into a vector
445 std::vector<unsigned int> runNumbers = JCT::utils::vectorize<unsigned int>(runNumbersString,", ");
446 if (runNumbers.empty())
447 {
448 ATH_MSG_ERROR("RunNumbers field specified, but value is empty or not unsigned ints for hot cells file: " << m_hotCellsFile);
449 return StatusCode::FAILURE;
450 }
451
452 // Loop over the run numbers
453 for (unsigned int run : runNumbers){
454 std::vector<std::unique_ptr<JCT::HotCell>>& cellVec = m_hotCellsMap[run];
455
456 // The number of hot cells should be below 100 for a given run...
457 for (size_t iCell = 0; iCell < 100; ++iCell)
458 {
459 const TString baseName = Form("Run%u.Cell%zu.",run,iCell);
460
461 // Read the cell info
462 const int layer = readCells.GetValue(baseName+"Layer", -1 );
463 const float minEta = readCells.GetValue(baseName+"EtaMin",-10.);
464 const float maxEta = readCells.GetValue(baseName+"EtaMax", 10.);
465 const float minPhi = readCells.GetValue(baseName+"PhiMin",-10.);
466 const float maxPhi = readCells.GetValue(baseName+"PhiMax", 10.);
467
468 // Check that the input makes sense
469 if (layer < 0 && minEta < -5 && maxEta > 5 && minPhi < -5 && maxPhi > 5)
470 continue;
471 if (layer < 0 || minEta < -5 || maxEta > 5 || minPhi < -5 || maxPhi > 5)
472 {
473 ATH_MSG_ERROR("Partially specified cell - please check the file: " << m_hotCellsFile);
474 ATH_MSG_ERROR(Form("Got Layer=%d, EtaMin=%f, EtaMax=%f, PhiMin=%f, PhiMax=%f",layer,minEta,maxEta,minPhi,maxPhi));
475 return StatusCode::FAILURE;
476 }
477 cellVec.emplace_back(std::make_unique<JCT::HotCell>(layer,minEta,maxEta,minPhi,maxPhi));
478 }
479
480 // Ensure we found the expected run
481 if (cellVec.empty())
482 {
483 ATH_MSG_ERROR("Specified that Run# " << run << " contains hot cells, but did not find any corresponding cells in the file: " << m_hotCellsFile);
484 return StatusCode::FAILURE;
485 }
486 }
487
488 // Done
489 return StatusCode::SUCCESS;
490}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Handle class for reading a decoration on an object.
Handle class for adding a decoration to an object.
const float m_etaMin
SG::ConstAccessor< std::vector< float > > m_ePerSamp
const float m_etaMax
const float m_phiMax
const float m_phiMin
bool jetAffectedByHotCell(const xAOD::Jet &jet) const
HotCell(const int layer, const float etaMin, const float etaMax, const float phiMin, const float phiMax)
asg::AcceptInfo m_accept
asg::AcceptData accept(const int isJetClean, const int fmaxIndex) const
The DFCommonJets decoration accept method.
static CleaningLevel getCutLevel(const std::string &)
Helpers for cut names.
JetCleaningTool(const std::string &name="JetCleaningTool")
Standard constructor.
Gaudi::Property< std::string > m_jetContainerName
SG::ConstAccessor< char > m_acc_jetClean
void missingVariable(const std::string &varName) const
Gaudi::Property< std::string > m_cutName
Name of the cut.
Gaudi::Property< bool > m_useDecorations
std::unordered_map< unsigned int, std::vector< std::unique_ptr< JCT::HotCell > > > m_hotCellsMap
bool containsHotCells(const xAOD::Jet &jet, const unsigned int runNumber) const
Hot cell checks.
virtual StatusCode decorate(const xAOD::JetContainer &jets) const override
Decorate a jet collection without otherwise modifying it.
static std::string getCutName(const CleaningLevel)
Gaudi::Property< std::string > m_hotCellsFile
Hot cells caching.
CleaningLevel m_cutLevel
CleaningLevel
Levels of cut.
Gaudi::Property< bool > m_doUgly
SG::WriteDecorHandleKey< xAOD::JetContainer > m_jetCleanKey
Gaudi::Property< bool > m_useLooseDecorForTightCut
virtual ~JetCleaningTool()
Standard destructor.
StatusCode readHotCells()
Hot cells reading helper.
virtual StatusCode initialize() override
Initialize method.
Helper class to provide constant type-safe access to aux data.
Handle class for adding a decoration to an object.
void setCutResult(const std::string &cutName, bool cutResult)
Set the result of a cut, based on the cut name (safer).
Definition AcceptData.h:134
bool getCutResult(const std::string &cutName) const
Get the result of a cut, based on the cut name (safer).
Definition AcceptData.h:98
Class mimicking the AthMessaging class from the offline software.
AsgMessaging(const std::string &name)
Constructor with a name.
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58
bool vectorize(const TString &str, const TString &sep, std::vector< T > &result)
Jet_v1 Jet
Definition of the current "jet version".
JetContainer_v1 JetContainer
Definition of the current "jet container version".
int run(int argc, char *argv[])