ATLAS Offline Software
Loading...
Searching...
No Matches
JetTileCorrectionTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// ROOT include(s):
6
7// EDM include(s):
9
10// Local include(s):
12
13#ifndef XAOD_STANDALONE // For now metadata is Athena-only
15#endif
16
23
24#include <iostream>
25#include <fstream>
26#include <sstream>
27#include "stdlib.h"
28#include <cmath>
29
30const double width=2.*M_PI/64.;
31
32//Validity ranges
33const float iov_pt_min = 20000.; //pt_min > 20GeV
34const float iov_pt_max = 3000000.; //pt_max > 3TeV
35const float iov_aeta_max = 2.8; //|eta| < 2.8
36
37namespace CP {
38
39 using namespace JTC;
40
41 static const SG::Decorator<unsigned int> dec_status("TileStatus");
42 static const SG::Decorator<float> dec_ptraw("Ptraw");
43
44
45 JetTileCorrectionTool :: JetTileCorrectionTool( const std::string& name )
47 m_RJET(0),
49 {
50 declareProperty("CorrectionFileName", m_rootFileName="JetTileCorrection/JetTile_pFile_010216.root", "Parametrization file");
51 declareProperty("MaskedRegionsMap", m_bd_dead_mapFile="JetTileCorrection/Tile_maskedDB_Run2.conf", "Masked regions DB file");
52 declareProperty("UserMaskedRegions", m_v_user_dead, "List of (ad-hoc) TileCal dead modules"); //array of extra dead modules
53 }
54
55 JetTileCorrectionTool :: ~JetTileCorrectionTool()
56 {
57 if ( m_rootFile ){
58 delete m_rootFile;
59 m_rootFile = nullptr;
60 }
61 }
62
63 StatusCode JetTileCorrectionTool :: initialize() {
64
65 // Greet the user:
66 ATH_MSG_INFO( "Initialising..." );
67
68#ifdef XAOD_STANDALONE
69 // Retrieve the event information (check if MC or data)
70 const xAOD::EventInfo* ei(nullptr);
71 ATH_CHECK( evtStore()->retrieve( ei, "EventInfo" ) );
72
74#else
75 // Retrieve the metadata (check if MC or data)
76 std::string projectName = "";
77 ATH_CHECK( AthAnalysisHelper::retrieveMetadata("/TagInfo", "project_name", projectName, inputMetaStore() ) );
78 if ( projectName == "IS_SIMULATION" ) m_isMC = true;
79 else if (projectName.compare(0, 4, "data") == 0 ) m_isMC = false;
80 ATH_MSG_INFO("Set up JetTileCorrectionTool -- this is MC? " << m_isMC);
81#endif
82
83
84 //set RJET
85 setRJET( (float) RJET );
86
87 //Load user dead regions map
89
90 //Load dead regions in DB map
91 loadDeadDB();
92
93 // Load the ROOT file
94 std::string fname = PathResolverFindCalibFile(m_rootFileName);
95 m_rootFile = TFile::Open( fname.c_str(), "READ" );
96
97 if ( !m_rootFile ) {
98 ATH_MSG_ERROR("Parametrization file " << m_rootFileName << "not found!");
99 return StatusCode::FAILURE;
100 }
101
102 //Load parametrizations
103 TString sub="";
104
105 m_pars_LB={};
106 m_pars_EB={};
107
108
109 for(unsigned int ieta=0; ieta < Pix_eta ; ieta++){
110 for(unsigned int iphi=0; iphi < Pix_phi ; iphi++){
111 sub = Form("_%d_%d",ieta,iphi);
112
113 //flatten 2D-etaphi-to-1D
114 m_pars_LB[ieta+iphi*Pix_eta]=(TH1F*)m_rootFile->Get("param_LB_fit"+sub);
115 m_pars_EB[ieta+iphi*Pix_eta]=(TH1F*)m_rootFile->Get("param_EB_fit"+sub);
116 m_pars_LB[ieta+iphi*Pix_eta]->SetDirectory(nullptr);
117 m_pars_EB[ieta+iphi*Pix_eta]->SetDirectory(nullptr);
118 }
119 }
120
121 //***
122 m_core_sys_LB=(TH1F*)m_rootFile->Get("sys_LB");
123 m_core_sys_EB=(TH1F*)m_rootFile->Get("sys_EB");
124 //***
125
126 //book number of pt bins in parametrization for later
127 m_NbinsPt = m_pars_LB[0]->GetNbinsX();
128
129
130 // set up for default running without systematics
132 ATH_MSG_ERROR("Loading the central value systematic set failed.");
133 return StatusCode::FAILURE;
134 }
135 // Add the affecting systematics to the global registry
137 if (registry.registerSystematics(*this) != StatusCode::SUCCESS){
138 ATH_MSG_ERROR("Unable to register the systematics");
139 return StatusCode::FAILURE;
140 }
141
142 ATH_MSG_DEBUG("Successfully initialized! ");
143
144 // Return gracefully:
145 return StatusCode::SUCCESS;
146 }
147
148
149 CorrectionCode JetTileCorrectionTool :: applyCorrection( xAOD::Jet& jet ){
150
151 //init decorations
152 dec_ptraw(jet) = jet.pt();
153 dec_status(jet) = (unsigned int) TS::UNKNOWN;
154
155 //check validity range of the correction
156 if( std::abs(jet.eta()) > iov_aeta_max ) return CorrectionCode::OutOfValidityRange;
157 if( jet.pt() < iov_pt_min || jet.pt() > iov_pt_max ) return CorrectionCode::OutOfValidityRange;
158
159 JTC::TS status = TS::GOOD;
160 if( loadAllModules(jet, status) != StatusCode::SUCCESS )
162
163 dec_status(jet) = (unsigned int) status; //save status decoration
164
165 // Nothing to do, return gracefully:
166 if(status == TS::GOOD) return CorrectionCode::Ok;
167
168 //get relative-pt factor to correct for
169 std::vector<float> cfactors = getCorrections(jet);
170
171 //if not correction need it, just leave
172 if(cfactors.empty()) return CorrectionCode::Ok;
173
174 // Redefine the jet 4vector by scaling both pt and mass
175 float newPt = jet.pt();
176 float newM = jet.m();
177
178 for(auto cf : cfactors){
179
180 newPt /= cf;
181 newM /= cf;
182 }
183
184 // Set the new jet 4vector
185 xAOD::JetFourMom_t newp4;
186 newp4.SetCoordinates(newPt, jet.eta(), jet.phi(), newM);
187
188 jet.setJetP4( newp4 );
189
190 ATH_MSG_DEBUG("JetTileCorrection applied.");
191
192 // Return gracefully:
193 return CorrectionCode::Ok;
194 }
195
196 CorrectionCode JetTileCorrectionTool :: correctedCopy( const xAOD::Jet& input, xAOD::Jet*& output) {
197
198 // A sanity check:
199 if( output ) {
200 ATH_MSG_WARNING( "Non-null pointer received. "
201 "There's a possible memory leak!" );
202 }
203
204 ATH_MSG_DEBUG("making the copy");
205
206 // Create the copy:
207 std::unique_ptr< xAOD::Jet > newobj( new xAOD::Jet( input ) );
208
209 // Apply the correction to it:
210 const CorrectionCode result = this->applyCorrection( *newobj );
212 ATH_MSG_ERROR( "Failed to apply correction to jet" );
213 } else {
214 output = newobj.release();
215 }
216
217 // Return the value from applyCorrection:
218 return result;
219 }
220
221 bool JetTileCorrectionTool :: isAffectedBySystematic( const SystematicVariation& systematic ) const
222 {
224 return sys.find (systematic) != sys.end ();
225 }
226
227
228 SystematicSet JetTileCorrectionTool :: affectingSystematics() const
229 {
231
232 result.insert(CP::SystematicVariation("JET_TILECORR_Uncertainty", 1));
233 result.insert(CP::SystematicVariation("JET_TILECORR_Uncertainty", -1));
234
235 return result;
236 }
237
238 SystematicSet JetTileCorrectionTool :: recommendedSystematics() const
239 {
240 return affectingSystematics();
241 }
242
243
244 StatusCode JetTileCorrectionTool :: applySystematicVariation ( const SystematicSet& systConfig ) {
245
246 // First, check if we already know this systematic configuration
247 auto itr = m_systFilter.find(systConfig);
248
249 // If it's a new input set, we need to filter it
250 if( itr == m_systFilter.end() ){
251
252 // New systematic. We need to parse it.
253 static const CP::SystematicSet affectingSys = affectingSystematics();
254 CP::SystematicSet filteredSys;
255 if (!CP::SystematicSet::filterForAffectingSystematics(systConfig, affectingSys, filteredSys)){
256 ATH_MSG_ERROR("Unsupported combination of systematics passed to the tool!");
257 return StatusCode::FAILURE;
258 }
259
260 // Insert filtered set into the map
261 itr = m_systFilter.insert(std::make_pair(systConfig, filteredSys)).first;
262 }
263
264 CP::SystematicSet& mySysConf = itr->second;
265 m_appliedSystematics = &mySysConf;
266 return StatusCode::SUCCESS;
267 }
268
269
270 JTC::TS JetTileCorrectionTool :: overlap(const xAOD::Jet& j, const Hole& region){
271
272 double phisize = region.phi2-region.phi1;
273 double phicenter = (region.phi1+region.phi2)/2.;
274 float jet_eta = j.jetP4(xAOD::JetConstitScaleMomentum).eta();
275 float jet_phi = j.jetP4(xAOD::JetConstitScaleMomentum).phi();
276
277 // GOOD?
278 if( region.eta2 < (jet_eta-m_RJET) ) return TS::GOOD;
279 if( region.eta1 > (jet_eta+m_RJET) ) return TS::GOOD;
280 if( std::abs(TVector2::Phi_mpi_pi(jet_phi-phicenter)) > m_RJET+phisize/2.) return TS::GOOD;
281
282 // CORE-BAD?
283 if( inHole(jet_eta, jet_phi, region) ) return TS::CORE;
284
285 // CORE-EDGE? (elsewhere)
286 return TS::EDGE;
287 }
288
289 void JetTileCorrectionTool :: loadDeadUser(){
290
291 m_user_dead_LB = {};
292 m_user_dead_EB = {};
293
294 //Format accepted for module input: "LBA 10" == "0 9"
295 for(auto& r : m_v_user_dead){
296
297 std::stringstream ss(r);
298 std::string s;
299
300 std::vector<std::string> tokens;
301 while (getline(ss, s, ' ')) {
302 tokens.push_back(s);
303 }
304
305 if(tokens.size() < 2){
306 ATH_MSG_ERROR("Part-Module pair " << r << " not known! Please use \"PART MOD\" format ");
307 continue;
308 }
309
310 int i_part = 0;
311 int i_mod = std::atoi((tokens.at(1)).c_str());
312 if (tokens.at(0).find('B') != std::string::npos)
313 i_mod -= 1; // substract 1 from second coor if given in format "LBA 4"
314
315 if(tokens.at(0)=="LBA" || tokens.at(0)== "0"){
316 i_part = 0;
317 }else if (tokens.at(0)=="LBC" || tokens.at(0)== "1"){
318 i_part = 1;
319 }else if (tokens.at(0)=="EBA" || tokens.at(0)== "2"){
320 i_part = 2;
321 }else if (tokens.at(0)=="EBC" || tokens.at(0)== "3"){
322 i_part = 3;
323 }else{
324 ATH_MSG_ERROR("Part-Module pair " << r << " not known! Please use \"PART MOD\" format ");
325 continue;
326 }
327
328 //Classify into long and extended barrel modules
329 Hole rdead = partModToHole(i_part, i_mod);
330
331 if( i_part < 2 )
332 m_user_dead_LB[r] = rdead; //LBA, LBC
333 else
334 m_user_dead_EB[r] = rdead; //EBA, EBC
335
336 }
337 }
338
339 //Load dead regions from DB file
340 void JetTileCorrectionTool :: loadDeadDB(){
341
342 m_db_dead_LB = {};
343 m_db_dead_EB = {};
344
345 std::vector<Hole> dbholes={};
346
347 //Simulation
348 if(m_isMC){
349
350 //-- MC15c : no dead modules
351 //COOLOFL_TILE/OFLP200 /TILE/OFL02/STATUS/ADC TileOfl02StatusAdc-IOVDEP-05
352 //COOLOFL_TILE/OFLP200 /TILE/OFL02/NOISE/CELL TileOfl02NoiseCell-OF2-07
353 //
354
355 //TODO: add RunNumber check and add these below if running on MC15b //M.T.
356
357 //-- MC15b : two dead modules
358 //COOLOFL_TILE/OFLP200 /TILE/OFL02/STATUS/ADC TileOfl02StatusAdc-IOVDEP-05
359 //COOLOFL_TILE/OFLP200 /TILE/OFL02/NOISE/CELL TileOfl02NoiseCell-OF2-07
360 //
361 // //LBA10
362 // Hole rdead = partModToHole(0, 9);
363 // rdead.iov = make_pair(0,1000000); //no range
364 // dbholes.push_back( rdead );
365
366 // //EBC21
367 // rdead = partModToHole(3, 20);
368 // rdead.iov = make_pair(0,1000000); //no range
369 // dbholes.push_back( rdead );
370
371 }
372 else{ //DATA
373
374 //read map file
375 std::string mapFilename = PathResolverFindCalibFile(m_bd_dead_mapFile);
376 std::ifstream mapFile;
377 mapFile.open(mapFilename);
378
379 std::string line;
380 while (std::getline(mapFile, line)){
381
382 if(line[0]=='#') continue;
383
384 std::istringstream iss(line);
385 int part,mod, irun, erun;
386 std::string modname;
387
388 if (!(iss >> part >> mod >> irun >> erun >> modname)) { break; } // error
389
390 Hole rdead = partModToHole(part, mod);
391 rdead.iov = std::make_pair(irun,erun);
392 dbholes.push_back( rdead );
393
394 }
395 }
396
397 int dbh=1;
398 for(const auto& h : dbholes){
399 if(std::abs(h.eta1)>1 || std::abs(h.eta2)>1){
400 m_db_dead_EB[Form("DB%d",dbh)] = h;
401 }
402 else{
403 m_db_dead_LB[Form("DB%d",dbh)] = h;
404 }
405
406 dbh++;
407
408 // ATH_MSG_INFO("Adding DB dead module at (eta1,phi1)=(" << h.eta1 << "," << h.phi1 << ")");
409 ATH_MSG_DEBUG("Adding DB dead module at (eta1,phi1)=(" << h.eta1 << "," << h.phi1 << ")");
410 }
411 }
412
413
414 IPair JetTileCorrectionTool :: getModulePosition(const xAOD::Jet& jet, const Hole& mod){
415
416 float eta = (mod.eta1+mod.eta2)/2.;
417 float phi = (mod.phi1+mod.phi2)/2.;
418
419 float jet_eta = jet.jetP4(xAOD::JetConstitScaleMomentum).eta();
420 float jet_phi = jet.jetP4(xAOD::JetConstitScaleMomentum).phi();
421
422 float eta_dist = jet_eta-eta;
423 float phi_dist = TVector2::Phi_mpi_pi((double)jet_phi - (double)phi);
424
425 int inphi;
426
427 if(std::abs(phi_dist)<0.05){ inphi=0; }
428 else if(std::abs(phi_dist)<0.1){ inphi=1; }
429 else if(std::abs(phi_dist)<0.2){ inphi=2; }
430 else if(std::abs(phi_dist)<0.3){ inphi=3; }
431 else if(std::abs(phi_dist)<0.4){ inphi=4; }
432 else { inphi=5;}
433
434 float ieta = eta_dist/PIXWIDTH;
435 int ineta = (int)ieta;
436
437
438 //get parametrization from positive side always!
439 if(mod.eta1 < -0.1) ineta = -ineta;
440
441 return std::make_pair(ineta+8, inphi);
442 }
443
444 void JetTileCorrectionTool :: loadModulesFromMap(const xAOD::Jet& jet, JTC::TS &status, const std::map<std::string,Hole>& hmap, PART part, TYPE type){
445
446 float cfactor(0.);
447 IPair mpos;
448 JTC::TS cstatus = TS::GOOD;
449
450 for (const auto& mod : hmap){
451
452 //check IOV for DB defined modules
453 if(!m_isMC && type == TYPE::DB && !inIOV(mod.second, m_current_run)) continue;
454
455 cstatus = overlap(jet,mod.second);
456 if(cstatus != TS::GOOD){
457
458 //get eta-phi position (relative to jet axis)
459 mpos = getModulePosition(jet,mod.second);
460
461 //get correction (1D-array)
462 if(part == PART::LB)
463 cfactor = m_pars_LB[mpos.first + Pix_eta * mpos.second]->GetBinContent(getPtBin(jet.pt()));
464 else
465 cfactor = m_pars_EB[mpos.first + Pix_eta * mpos.second]->GetBinContent(getPtBin(jet.pt()));
466
467 //overwrite global status with worst case
468 if(cstatus > status) status = cstatus;
469
470 //create new d-region and add it to the list
471 m_position_masked.emplace_back( mod.second,mpos,part,cfactor,cstatus,type );
472
473 }
474 }
475
476 }
477
478
479 JTC::TS JetTileCorrectionTool :: getTileStatus(const xAOD::Jet& jet){
480 JTC::TS status = TS::GOOD;
481 if( loadAllModules(jet, status) != StatusCode::SUCCESS ){
482 ATH_MSG_ERROR( "Something went wrong while loading/checking the modules!");
483 return TS::UNKNOWN;
484 }
485 return status;
486 }
487
488
489 StatusCode JetTileCorrectionTool :: addTileStatus(const xAOD::Jet& jet){
490 JTC::TS status = getTileStatus(jet);
491 dec_status(jet) = (unsigned int) status; //save status decoration
492 return StatusCode::SUCCESS;
493 }
494
495 void JetTileCorrectionTool :: setRJET(float r){
496 m_RJET = r;
497 }
498
499
500 StatusCode JetTileCorrectionTool :: loadAllModules(const xAOD::Jet& jet, JTC::TS &status){
501
502 m_position_masked.clear();
503
504 //Read dead regions DB if we change run/LB
505 if(!m_isMC){
506
507 const xAOD::EventInfo* ei = nullptr;
508 if( evtStore()->retrieve( ei, "EventInfo" ).isFailure() ) {
509 ATH_MSG_WARNING( "No EventInfo object could be retrieved" );
510 }
511
512 if (!ei) {
513 ATH_MSG_ERROR( "Cannot retrieve the EventInfo" );
514 return StatusCode::FAILURE;
515 }
516
517 m_current_run = ei->runNumber();
518
519 }
520
521 //load DB-defined modules in LB
523 //load DB-defined modules in EB
525
526 //load user-defined modules in LB
528 //load user-defined modules in EB
530
531 //sort modules by correction size (in decreasing order)
532 std::sort(m_position_masked.begin(), m_position_masked.end(), std::greater<Region>());
533
534 return StatusCode::SUCCESS;
535 }
536
537 std::vector<float> JetTileCorrectionTool :: getCorrections(const xAOD::Jet& jet){
538
539 float ptlast = jet.pt();
540 int ptbin = getPtBin(ptlast);
541
542 //***
543 int ptbin_sys=m_core_sys_LB->FindBin(ptlast);
544 //***
545
546 std::vector<float> corrections = {};
547
548 if (ptbin < 0) return corrections; // no correction below 20GeV
549
550 //***
551 if (ptlast < 40000) return corrections;
552 //***
553
554 float clast = 1.;
555 float sigma = 0.;
556
557 if(m_appliedSystematics!=nullptr){ //nomimal value
558 sigma = m_appliedSystematics->getParameterByBaseName("JET_TILECORR_Uncertainty");
559 }
560
561 for(const auto& region : m_position_masked){
562
563 ptlast *= clast;
564 ptbin = getPtBin( ptlast );
565
566 if(region.part == PART::LB) {
567
568 //***
569 if (region.status == TS::CORE) clast = m_pars_LB[region.ep.first + Pix_eta * region.ep.second]->GetBinContent(ptbin) + sigma * m_core_sys_LB->GetBinContent(ptbin_sys); //systematics are bigger if we are in the core
570 else clast = m_pars_LB[region.ep.first + Pix_eta * region.ep.second]->GetBinContent(ptbin) + sigma * m_pars_LB[region.ep.first + Pix_eta * region.ep.second]->GetBinError(ptbin);
571 //***
572
573 }else{
574
575 //***
576 if (region.status == TS::CORE) clast = m_pars_EB[region.ep.first + Pix_eta * region.ep.second]->GetBinContent(ptbin) + sigma * m_core_sys_EB->GetBinContent(ptbin_sys); //systematics are bigger if we are in the core
577 else clast = m_pars_EB[region.ep.first + Pix_eta * region.ep.second]->GetBinContent(ptbin) + sigma * m_pars_EB[region.ep.first + Pix_eta * region.ep.second]->GetBinError(ptbin);
578 //***
579
580
581 }
582 corrections.push_back(clast);
583 }
584
585 return corrections;
586 }
587
588
589 int JetTileCorrectionTool :: getPtBin(float pt){
590
591 if(pt < 20000.) return -1;
592
593 int ptbin = m_pars_LB[0]->FindBin(pt*0.001); // no correction below 20GeV
594
595 //adjust under/overflows
596 if (ptbin < 1) ptbin = 1;
597 if (ptbin > m_NbinsPt) ptbin = m_NbinsPt;
598
599 return ptbin;
600 }
601
602 bool JetTileCorrectionTool :: inIOV(const Hole& region, int run){
603 if( region.iov.first>=0 && (run < region.iov.first) ) return false;
604 if( region.iov.second>=0 && (run > region.iov.second) ) return false;
605 return true;
606 }
607
608 bool JetTileCorrectionTool :: inHole(const xAOD::Jet& jet, const Hole& rdead){
609 return inHole(jet.eta(), jet.phi(), rdead);
610 }
611
612
613 bool JetTileCorrectionTool :: inHole(float eta, float phi, const Hole& rdead){
614
615 if(rdead.eta1==rdead.eta2 || rdead.phi1==rdead.phi2) return false;
616 if((eta > rdead.eta1) && (eta < rdead.eta2) && (phi > rdead.phi1) && (phi < rdead.phi2)) return true;
617
618 return false;
619 }
620
622
623 Hole region;
624 region.iov = std::make_pair(0,1000000); //a dummy full IOV is set by default
625
626 switch(part){
627 case 0: //LBA
628 region.eta1=0;
629 region.eta2=0.9;
630 break;
631 case 1: //LBC
632 region.eta1=-0.9;
633 region.eta2=0;
634 break;
635 case 2: //EBA
636 region.eta1=0.8;
637 region.eta2=1.7;
638 break;
639 case 3: //EBC
640 region.eta1=-1.7;
641 region.eta2=-0.8;
642 break;
643 default:
644 std::cout<<"Bad partition value passed!\n";
645 region.eta1=-999;
646 region.eta2=-999;
647 region.phi1=-999;
648 region.phi2=-999;
649 return region;
650 }
651 if(mod<32){
652 region.phi1=((double)mod)*width;
653 region.phi2=region.phi1+width;
654 }
655 else if(mod<64){
656 region.phi1=((double)mod)*width-2.*M_PI;
657 region.phi2=region.phi1+width;
658 }
659 else{
660 std::cout<<"Bad module value passed!\n";
661 region.eta1=-999;
662 region.eta2=-999;
663 region.phi1=-999;
664 region.phi2=-999;
665 return region;
666 }
667 return region;
668 }
669
670
671} // namespace CP
#define M_PI
Scalar eta() const
pseudorapidity method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
Helper class to provide type-safe access to aux data.
std::pair< int, int > IPair
const float iov_aeta_max
const float iov_pt_max
const float iov_pt_min
#define Pix_phi
#define Pix_eta
#define RJET
#define PIXWIDTH
static Double_t ss
std::string PathResolverFindCalibFile(const std::string &logical_file_name)
const double width
static std::string retrieveMetadata(const std::string &folder, const std::string &key, const ServiceHandle< StoreGateSvc > &inputMetaStore)
method that always returns as a string you can use from, e.g, pyROOT with evt = ROOT....
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
ServiceHandle< StoreGateSvc > & evtStore()
Header file for AthHistogramAlgorithm.
Return value from object correction CP tools.
@ Error
Some error happened during the object correction.
@ OutOfValidityRange
Input object is out of validity range.
@ Ok
The correction was done successfully.
TH1F * m_core_sys_LB
//needs to be transient!
std::map< std::string, JTC::Hole > m_db_dead_LB
std::map< std::string, JTC::Hole > m_user_dead_LB
IPair getModulePosition(const xAOD::Jet &jet, const JTC::Hole &module)
bool inIOV(const JTC::Hole &region, int run)
bool inHole(const xAOD::Jet &j, const JTC::Hole &rdead)
std::map< int, TH1F * > m_pars_EB
virtual StatusCode applySystematicVariation(const SystematicSet &systConfig)
effects: configure this tool for the given list of systematic variations.
virtual CorrectionCode applyCorrection(xAOD::Jet &object)
Apply the correction on a modifyable object.
std::vector< JTC::Region > m_position_masked
std::vector< std::string > m_v_user_dead
std::unordered_map< CP::SystematicSet, CP::SystematicSet > m_systFilter
Systematics filter map.
std::map< int, TH1F * > m_pars_LB
JTC::Hole partModToHole(int part, int mod)
std::map< std::string, JTC::Hole > m_user_dead_EB
StatusCode loadAllModules(const xAOD::Jet &jet, JTC::TS &status)
CP::SystematicSet * m_appliedSystematics
Currently applied systematics.
JTC::TS overlap(const xAOD::Jet &j, const JTC::Hole &region)
void loadModulesFromMap(const xAOD::Jet &jet, JTC::TS &status, const std::map< std::string, JTC::Hole > &hmap, JTC::PART part=JTC::PART::LB, JTC::TYPE type=JTC::TYPE::DB)
virtual SystematicSet affectingSystematics() const
returns: the list of all systematics this tool can be affected by
JTC::TS getTileStatus(const xAOD::Jet &jet)
std::map< std::string, JTC::Hole > m_db_dead_EB
This module implements the central registry for handling systematic uncertainties with CP tools.
static SystematicRegistry & getInstance()
Get the singleton instance of the registry for the curren thread.
StatusCode registerSystematics(const IReentrantSystematicsTool &tool)
effects: register all the systematics from the tool
Class to wrap a set of SystematicVariations.
static StatusCode filterForAffectingSystematics(const SystematicSet &systConfig, const SystematicSet &affectingSystematics, SystematicSet &filteredSystematics)
description: filter the systematics for the affected systematics returns: success guarantee: strong f...
Helper class to provide type-safe access to aux data.
Definition Decorator.h:59
AsgMetadataTool(const std::string &name)
Normal ASG tool constructor with a name.
MetaStorePtr_t inputMetaStore() const
Accessor for the input metadata store.
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.
JetFourMom_t jetP4() const
The full 4-momentum of the particle : internal jet type.
Definition Jet_v1.cxx:76
int r
Definition globals.cxx:22
Select isolated Photons, Electrons and Muons.
static const SG::Decorator< float > dec_ptraw("Ptraw")
static const SG::Decorator< unsigned int > dec_status("TileStatus")
Definition run.py:1
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
Jet_v1 Jet
Definition of the current "jet version".
setRcore setEtHad setFside pt
EventInfo_v1 EventInfo
Definition of the latest event info version.
setRawEt setRawPhi int
@ JetConstitScaleMomentum
Definition JetTypes.h:29
ROOT::Math::LorentzVector< ROOT::Math::PtEtaPhiM4D< double > > JetFourMom_t
Base 4 Momentum type for Jet.
Definition JetTypes.h:17
double getCorrections(TH2F *corrections, double eta, double phi, double Et, double charge)
std::pair< int, int > iov