ATLAS Offline Software
Loading...
Searching...
No Matches
CommonEfficiencyTool.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// Framework include(s):
7
8// local include(s)
11
12// ROOT include(s)
13#include "TF1.h"
14#include "TH1.h"
15#include "TH2.h"
16#include "TROOT.h"
17#include "TClass.h"
18
19using namespace TauAnalysisTools;
20
21/*
22 This tool acts as a common tool to apply efficiency scale factors and
23 uncertainties. By default, only nominal scale factors without systematic
24 variations are applied. Unavailable systematic variations are ignored, meaning
25 that the tool only returns the nominal value. In case the one available
26 systematic is requested, the smeared scale factor is computed as:
27 - sf = sf_nominal +/- n * uncertainty
28
29 where n is in general 1 (representing a 1 sigma smearing), but can be any
30 arbitrary value. In case multiple systematic variations are passed they are
31 added in quadrature. Note that it's currently only supported if all are up or
32 down systematics.
33
34 The tool reads in root files including TH2 histograms which need to fullfil a
35 predefined structure:
36
37 scale factors:
38 - sf_<workingpoint>_<prongness>p
39 uncertainties:
40 - <NP>_<up/down>_<workingpoint>_<prongness>p (for asymmetric uncertainties)
41 - <NP>_<workingpoint>_<prongness>p (for symmetric uncertainties)
42
43 where the <workingpoint> (e.g. loose/medium/tight) fields may be
44 optional. <prongness> represents either 1 or 3, whereas 3 is currently used
45 for multiprong in general. The <NP> fields are names for the type of nuisance
46 parameter (e.g. STAT or SYST), note the tool decides whethe the NP is a
47 recommended or only an available systematic based on the first character:
48 - uppercase -> recommended
49 - lowercase -> available
50 This magic happens here:
51 - CommonEfficiencyTool::generateSystematicSets()
52
53 In addition the root input file can also contain objects of type TF1 that can
54 be used to provide kind of unbinned scale factors or systematics. The major
55 usecase for now is the high-pt uncertainty for the tau ID and tau
56 reconstruction.
57
58 The files may also include TNamed objects which is used to define how x and
59 y-axes should be treated. By default the x-axis is given in units of tau-pT in
60 GeV and the y-axis is given as tau-eta. If there is for example a TNamed
61 object with name "Yaxis" and title "|eta|" the y-axis is treated in units of
62 absolute tau eta. All this is done in:
63 - void CommonEfficiencyTool::ReadInputs(TFile* fFile)
64
65*/
66
67//______________________________________________________________________________
69 : asg::AsgTool( sName )
70 , m_mSF(nullptr)
71 , m_sSystematicSet(nullptr)
74 , m_sSFHistName("sf")
75 , m_bNoMultiprong(false)
77 , m_bSFIsAvailable(false)
79{
80}
81
82/*
83 need to clear the map of histograms cause we have the ownership, not ROOT
84*/
86{
87 if (m_mSF)
88 for (auto mEntry : *m_mSF)
89 delete std::get<0>(mEntry.second);
90}
91
92/*
93 - Find the root files with scale factor inputs on cvmfs using PathResolver
94 (more info here:
95 https://twiki.cern.ch/twiki/bin/viewauth/AtlasComputing/PathResolver)
96 - Call further functions to process and define NP strings and so on
97 - Configure to provide nominal scale factors by default
98*/
100{
101 ATH_MSG_INFO( "Initializing CommonEfficiencyTool" );
102 // only read in histograms once
103 if (m_mSF==nullptr)
104 {
105 std::string sInputFilePath = PathResolverFindCalibFile(m_sInputFilePath);
106
107 m_mSF = std::make_unique< tSFMAP >();
108 std::unique_ptr< TFile > fSF( TFile::Open( sInputFilePath.c_str(), "READ" ) );
109 if(!fSF)
110 {
111 ATH_MSG_FATAL("Could not open file " << sInputFilePath.c_str());
112 return StatusCode::FAILURE;
113 }
114 ReadInputs(*fSF);
115 fSF->Close();
116 }
117
118 // needed later on in generateSystematicSets(), maybe move it there
119 std::vector<std::string> vInputFilePath;
120 split(m_sInputFilePath,'/',vInputFilePath);
121 m_sInputFileName = vInputFilePath.back();
122
124
125 if (!m_sWP.empty())
126 m_sSFHistName = "sf_"+m_sWP;
127
128 // load empty systematic variation by default
129 if (applySystematicVariation(CP::SystematicSet()) != StatusCode::SUCCESS )
130 return StatusCode::FAILURE;
131
132 return StatusCode::SUCCESS;
133}
134
135/*
136 Retrieve the scale factors and if requested the values for the NP's and add
137 this stuff in quadrature. Finally return sf_nom +/- n*uncertainty
138*/
139
140//______________________________________________________________________________
142 double& dEfficiencyScaleFactor, unsigned int /*iRunNumber*/)
143{
144 // check which true state is requested
146 {
147 dEfficiencyScaleFactor = 1.;
149 }
150
151 // check if 1 prong
152 if (m_bNoMultiprong && xTau.nTracks() != 1)
153 {
154 dEfficiencyScaleFactor = 1.;
156 }
157
158 // get decay mode or prong extension for histogram name
159 std::string sMode;
161 {
162 int iDecayMode = -1;
164 sMode = ConvertDecayModeToString(iDecayMode);
165 if (sMode.empty())
166 {
167 ATH_MSG_WARNING("Found tau with unknown decay mode. Skip efficiency correction.");
169 }
170 }
171 else
172 {
173 // skip taus which are not 1 or 3 prong
174 if( xTau.nTracks() != 1 && xTau.nTracks() != 3) {
175 dEfficiencyScaleFactor = 1.;
177 }
178
179 sMode = ConvertProngToString(xTau.nTracks());
180 }
181
182 std::string sHistName;
183 if(m_doTauTrig){
184 sHistName = "sf_all_"+m_sWP+sMode;
185 } else {
186 sHistName = m_sSFHistName + sMode;
187 }
188
189 // get standard scale factor
190 CP::CorrectionCode tmpCorrectionCode = getValue(sHistName,
191 xTau,
192 dEfficiencyScaleFactor);
193 // return correction code if histogram is not available
194 if (tmpCorrectionCode != CP::CorrectionCode::Ok)
195 return tmpCorrectionCode;
196
197 // skip further process if systematic set is empty
198 if (m_sSystematicSet->empty())
200
201 // get uncertainties summed in quadrature
202 double dTotalSystematic2 = 0.;
203 double dDirection = 0.;
204 for (auto syst : *m_sSystematicSet)
205 {
206 // check if systematic is available
207 auto it = m_mSystematicsHistNames.find(syst.basename());
208
209 // get uncertainty value
210 double dUncertaintySyst = 0.;
211
212 // needed for up/down decision
213 dDirection = syst.parameter();
214
215 // build up histogram name
216 sHistName = it->second;
217 if (dDirection>0.) sHistName+="_up";
218 else sHistName+="_down";
219
220 if(m_doTauTrig){ sHistName+="_all"; }
221
222 if (!m_sWP.empty()) sHistName+="_"+m_sWP;
223 sHistName += sMode;
224
225
226 // filter unwanted combinations
227 if( (sHistName.find("3P") != std::string::npos && sHistName.find("1p") != std::string::npos) ||
228 (sHistName.find("1P") != std::string::npos && sHistName.find("3p") != std::string::npos))
229 continue;
230
231 if( (sHistName.find("1520") != std::string::npos && sHistName.find("loose") != std::string::npos) ){
232 continue;
233 }
234
235 // get the uncertainty from the histogram
236 tmpCorrectionCode = getValue(sHistName,
237 xTau,
238 dUncertaintySyst);
239
240 // return correction code if histogram is not available
241 if (tmpCorrectionCode != CP::CorrectionCode::Ok)
242 return tmpCorrectionCode;
243
244 // scale uncertainty with direction, i.e. +/- n*sigma
245 dUncertaintySyst *= dDirection;
246
247 // square uncertainty and add to total uncertainty
248 dTotalSystematic2 += dUncertaintySyst * dUncertaintySyst;
249 }
250
251 // now use dDirection to use up/down uncertainty
252 dDirection = (dDirection > 0.) ? 1. : -1.;
253
254 // finally apply uncertainty (eff * ( 1 +/- \sum )
255 dEfficiencyScaleFactor *= 1. + dDirection * std::sqrt(dTotalSystematic2);
256
258}
259
260/*
261 Get scale factor from getEfficiencyScaleFactor and decorate it to the
262 tau. Note that this can only be done if the variable name is not already used,
263 e.g. if the variable was already decorated on a previous step (enured by the
264 m_bSFIsAvailableChecked check).
265
266 Technical note: cannot use `static SG::Decorator` as we will have
267 multiple instances of this tool with different decoration names.
268*/
269//______________________________________________________________________________
271 unsigned int iRunNumber)
272{
273 double dSf = 0.;
274
277 {
278 m_bSFIsAvailable = decor.isAvailable(xTau);
281 {
282 ATH_MSG_DEBUG(m_sVarName << " decoration is available on first tau processed, switched off applyEfficiencyScaleFactor for further taus.");
283 ATH_MSG_DEBUG("If an application of efficiency scale factors needs to be redone, please pass a shallow copy of the original tau.");
284 }
285 }
288
289 // retrieve scale factor
290 CP::CorrectionCode tmpCorrectionCode = getEfficiencyScaleFactor(xTau, dSf, iRunNumber);
291 // adding scale factor to tau as decoration
292 decor(xTau) = dSf;
293
294 return tmpCorrectionCode;
295}
296
297/*
298 standard check if a systematic is available
299*/
300//______________________________________________________________________________
302{
304 return sys.find(systematic) != sys.end();
305}
306
307/*
308 standard way to return systematics that are available (including recommended
309 systematics)
310*/
311//______________________________________________________________________________
316
317/*
318 standard way to return systematics that are recommended
319*/
320//______________________________________________________________________________
325
326/*
327 Configure the tool to use a systematic variation for further usage, until the
328 tool is reconfigured with this function. The passed systematic set is checked
329 for sanity:
330 - unsupported systematics are skipped
331 - only combinations of up or down supported systematics is allowed
332 - don't mix recommended systematics with other available systematics, cause
333 sometimes recommended are a quadratic sum of the other variations,
334 e.g. TOTAL=(SYST^2 + STAT^2)^0.5
335*/
336//______________________________________________________________________________
338{
339
340 // first check if we already know this systematic configuration
341 auto itSystematicSet = m_mSystematicSets.find(sSystematicSet);
342 if (itSystematicSet != m_mSystematicSets.end())
343 {
344 m_sSystematicSet = &itSystematicSet->first;
345 return StatusCode::SUCCESS;
346 }
347
348 // sanity checks if systematic set is supported
349 double dDirection = 0.;
350 CP::SystematicSet sSystematicSetAvailable;
351 for (auto sSyst : sSystematicSet)
352 {
353 // check if systematic is available
354 auto it = m_mSystematicsHistNames.find(sSyst.basename());
355 if (it == m_mSystematicsHistNames.end())
356 {
357 ATH_MSG_VERBOSE("unsupported systematic variation: "<< sSyst.basename()<<"; skipping this one");
358 continue;
359 }
360
361
362 if (sSyst.parameter() * dDirection < 0)
363 {
364 ATH_MSG_ERROR("unsupported set of systematic variations, you should either use only \"UP\" or only \"DOWN\" systematics in one set!");
365 ATH_MSG_ERROR("systematic set will not be applied");
366 return StatusCode::FAILURE;
367 }
368 dDirection = sSyst.parameter();
369
370 if ((m_sRecommendedSystematics.find(sSyst.basename()) != m_sRecommendedSystematics.end()) and sSystematicSet.size() > 1)
371 {
372 ATH_MSG_ERROR("unsupported set of systematic variations, you should not combine \"TAUS_{TRUE|FAKE}_EFF_*_TOTAL\" with other systematic variations!");
373 ATH_MSG_ERROR("systematic set will not be applied");
374 return StatusCode::FAILURE;
375 }
376
377 // finally add the systematic to the set of systematics to process
378 sSystematicSetAvailable.insert(sSyst);
379 }
380
381 // store this calibration for future use, and make it current
382 m_sSystematicSet = &m_mSystematicSets.insert(std::pair<CP::SystematicSet,std::string>(sSystematicSetAvailable, sSystematicSet.name())).first->first;
383
384 return StatusCode::SUCCESS;
385}
386
387//=================================PRIVATE-PART=================================
388std::string CommonEfficiencyTool::ConvertProngToString(const int fProngness) const
389{
390 return fProngness == 1 ? "_1p" : "_3p";
391}
392
393/*
394 decay mode converter
395*/
396//______________________________________________________________________________
397std::string CommonEfficiencyTool::ConvertDecayModeToString(const int iDecayMode) const
398{
399 switch(iDecayMode)
400 {
402 return "_r1p0n";
404 return "_r1p1n";
406 return "_r1pXn";
408 return "_r3p0n";
410 return "_r3pXn";
411 default:
412 return "";
413 }
414}
415
416/*
417 Read in a root file and store all objects to a map of this type:
418 std::map<std::string, tTupleObjectFunc > (see header) It's basically a map of
419 the histogram name and a function pointer based on the TObject type (TH1F,
420 TH1D, TF1). This is resolved in the function:
421 - CommonEfficiencyTool::addHistogramToSFMap
422 Further this function figures out the axis definition (see description on the
423 top)
424*/
425//______________________________________________________________________________
426void CommonEfficiencyTool::ReadInputs(const TFile& fFile)
427{
428 m_mSF->clear();
429
430 // initialize function pointer
431 m_fX = &finalTauPt;
432 m_fY = &finalTauEta;
433
434 TKey *kKey;
435 TIter itNext(fFile.GetListOfKeys());
436 while ((kKey = (TKey*)itNext()))
437 {
438 // parse file content for objects of type TNamed, check their title for
439 // known strings and reset funtion pointer
440 std::string sKeyName = kKey->GetName();
441 if (sKeyName == "Xaxis")
442 {
443 TNamed* tObj = (TNamed*)kKey->ReadObj();
444 std::string sTitle = tObj->GetTitle();
445 delete tObj;
446 if (sTitle == "TruthDecayMode")
447 {
449 ATH_MSG_DEBUG("using truth decay mode for x-axis");
450 }
451 if (sTitle == "truth visible pt")
452 {
454 ATH_MSG_DEBUG("using truth visible pT for x-axis");
455 }
456 if (sTitle == "|eta|")
457 {
459 ATH_MSG_DEBUG("using absolute tau eta for x-axis");
460 }
461
462 continue;
463 }
464 else if (sKeyName == "Yaxis")
465 {
466 TNamed* tObj = (TNamed*)kKey->ReadObj();
467 std::string sTitle = tObj->GetTitle();
468 delete tObj;
469 if (sTitle == "|eta|")
470 {
472 ATH_MSG_DEBUG("using absolute tau eta for y-axis");
473 }
474 else if (sTitle == "truth |eta|")
475 {
477 ATH_MSG_DEBUG("using absolute truth tau eta for y-axis");
478 }
479 continue;
480 }
481
482 std::vector<std::string> vSplitName = {};
483 split(sKeyName,'_',vSplitName);
484 if (vSplitName[0] == "sf")
485 {
486 addHistogramToSFMap(kKey, sKeyName);
487 }
488 else
489 {
490 // std::string sDirection = vSplitName[1];
491 if (sKeyName.find("_up_") != std::string::npos or sKeyName.find("_down_") != std::string::npos)
492 addHistogramToSFMap(kKey, sKeyName);
493 else
494 {
495 size_t iPos = sKeyName.find('_');
496 addHistogramToSFMap(kKey, sKeyName.substr(0,iPos)+"_up"+sKeyName.substr(iPos));
497 addHistogramToSFMap(kKey, sKeyName.substr(0,iPos)+"_down"+sKeyName.substr(iPos));
498 }
499 }
500 }
501 ATH_MSG_INFO("data loaded from " << fFile.GetName());
502}
503
504/*
505 Create the tuple objects for the map
506*/
507//______________________________________________________________________________
508void CommonEfficiencyTool::addHistogramToSFMap(TKey* kKey, const std::string& sKeyName)
509{
510 // handling for the 3 different input types TH1F/TH1D/TF1, function pointer
511 // handle the access methods for the final scale factor retrieval
512 TClass *cClass = gROOT->GetClass(kKey->GetClassName());
513 if (cClass->InheritsFrom("TH2"))
514 {
515 TH1* oObject = (TH1*)kKey->ReadObj();
516 oObject->SetDirectory(0);
517 (*m_mSF)[sKeyName] = tTupleObjectFunc(oObject,&getValueTH2);
518 ATH_MSG_DEBUG("added histogram with name "<<sKeyName);
519 }
520 else if (cClass->InheritsFrom("TH1"))
521 {
522 TH1* oObject = (TH1*)kKey->ReadObj();
523 oObject->SetDirectory(0);
524 (*m_mSF)[sKeyName] = tTupleObjectFunc(oObject,&getValueTH1);
525 ATH_MSG_DEBUG("added histogram with name "<<sKeyName);
526 }
527 else if (cClass->InheritsFrom("TF1"))
528 {
529 TObject* oObject = kKey->ReadObj();
530 (*m_mSF)[sKeyName] = tTupleObjectFunc(oObject,&getValueTF1);
531 ATH_MSG_DEBUG("added function with name "<<sKeyName);
532 }
533 else
534 {
535 ATH_MSG_DEBUG("ignored object with name "<<sKeyName);
536 }
537}
538
539/*
540 This function parses the names of the obejects from the input file and
541 generates the systematic sets and defines which ones are recommended or only
542 available. It also checks, based on the root file name, on which tau it needs
543 to be applied, e.g. only on reco taus coming from true taus or on those faked
544 by true electrons...
545
546 Examples:
547 filename: Reco_TrueHadTau_2016-ichep.root -> apply only to true taus
548 histname: sf_1p -> nominal 1p scale factor
549 histname: TOTAL_3p -> "total" 3p NP, recommended
550 histname: afii_1p -> "total" 3p NP, not recommended, but available
551*/
552//______________________________________________________________________________
554{
555 // creation of basic string for all NPs, e.g. "TAUS_TRUEHADTAU_EFF_RECO_"
556 std::vector<std::string> vSplitInputFilePath = {};
557 split(m_sInputFileName,'_',vSplitInputFilePath);
558 std::string sEfficiencyType = vSplitInputFilePath.at(0);
559 std::string sTruthType = vSplitInputFilePath.at(1);
560 std::transform(sEfficiencyType.begin(), sEfficiencyType.end(), sEfficiencyType.begin(), toupper);
561 std::transform(sTruthType.begin(), sTruthType.end(), sTruthType.begin(), toupper);
562 std::string sSystematicBaseString = "TAUS_"+sTruthType+"_EFF_"+sEfficiencyType+"_";
563
564 // set truth type to check for in truth matching
565 if (sTruthType=="TRUEHADTAU") m_eCheckTruth = TauAnalysisTools::TruthHadronicTau;
566 else if (sTruthType=="TRUEELECTRON") m_eCheckTruth = TauAnalysisTools::TruthElectron;
567 // 3p eVeto, still need this to be measurable in T&P
568 if (sEfficiencyType=="ELERNN" || sEfficiencyType=="ELEOLR") m_bNoMultiprong = true;
569
570 for (auto mSF : *m_mSF)
571 {
572 // parse for nuisance parameter in histogram name
573 std::vector<std::string> vSplitNP = {};
574 split(mSF.first,'_',vSplitNP);
575 std::string sNP = vSplitNP.at(0);
576 std::string sNPUppercase = vSplitNP.at(0);
577
578 // skip nominal scale factors
579 if (sNP == "sf") continue;
580
581 // skip if 3p histogram to avoid duplications (TODO: come up with a better solution)
582 //if (mSF.first.find("_3p") != std::string::npos) continue;
583
584 // test if NP starts with a capital letter indicating that this should be recommended
585 bool bIsRecommended = false;
586 if (isupper(sNP.at(0)) || isupper(sNP.at(1)))
587 bIsRecommended = true;
588
589 // make sNP uppercase and build final NP entry name
590 std::transform(sNPUppercase.begin(), sNPUppercase.end(), sNPUppercase.begin(), toupper);
591 std::string sSystematicString = sSystematicBaseString+sNPUppercase;
592
593 // add all found systematics to the AffectingSystematics
594 m_sAffectingSystematics.insert(CP::SystematicVariation (sSystematicString, 1));
595 m_sAffectingSystematics.insert(CP::SystematicVariation (sSystematicString, -1));
596 // only add found uppercase systematics to the RecommendedSystematics
597 if (bIsRecommended)
598 {
599 m_sRecommendedSystematics.insert(CP::SystematicVariation (sSystematicString, 1));
600 m_sRecommendedSystematics.insert(CP::SystematicVariation (sSystematicString, -1));
601 }
602
603 ATH_MSG_DEBUG("connected base name " << sNP << " with systematic " <<sSystematicString);
604 m_mSystematicsHistNames.insert({sSystematicString,sNP});
605 }
606}
607
608/*
609 return value from the tuple map object based on the pt/eta values (or the
610 corresponding value in case of configuration)
611*/
612//______________________________________________________________________________
614 const xAOD::TauJet& xTau,
615 double& dEfficiencyScaleFactor) const
616{
617 const tSFMAP& mSF = *m_mSF;
618 auto it = mSF.find (sHistName);
619 if (it == mSF.end())
620 {
621 ATH_MSG_ERROR("Object with name "<<sHistName<<" was not found in input file.");
622 ATH_MSG_DEBUG("Content of input file");
623 for (auto eEntry : mSF)
624 ATH_MSG_DEBUG(" Entry: "<<eEntry.first);
626 }
627
628 // get a tuple (TObject*,functionPointer) from the scale factor map
629 tTupleObjectFunc tTuple = it->second;
630
631 // get pt and eta (for x and y axis respectively)
632 double dPt = m_fX(xTau);
633 double dEta = m_fY(xTau);
634
635 double dVars[2] = {dPt, dEta};
636
637 // finally obtain efficiency scale factor from TH1F/TH1D/TF1, by calling the
638 // function pointer stored in the tuple from the scale factor map
639 return (std::get<1>(tTuple))(std::get<0>(tTuple), dEfficiencyScaleFactor, dVars);
640}
641
642/*
643 find the particular value in TH1 depending on pt (or the
644 corresponding value in case of configuration)
645 Note: In case values are outside of bin ranges, the closest bin value is used
646*/
647//______________________________________________________________________________
649 double& dEfficiencyScaleFactor, double dVars[])
650{
651 double dPt = dVars[0];
652
653 const TH1* hHist = dynamic_cast<const TH1*>(oObject);
654
655 if (!hHist)
656 {
657 // ATH_MSG_ERROR("Problem with casting TObject of type "<<oObject->ClassName()<<" to TH2F");
659 }
660
661 // protect values from underflow bins
662 dPt = std::max(dPt,hHist->GetXaxis()->GetXmin());
663 // protect values from overflow bins (times .999 to keep it inside last bin)
664 dPt = std::min(dPt,hHist->GetXaxis()->GetXmax() * .999);
665
666 // get bin from TH2 depending on x and y values; finally set the scale factor
667 int iBin = hHist->FindFixBin(dPt);
668 dEfficiencyScaleFactor = hHist->GetBinContent(iBin);
670}
671
672/*
673 find the particular value in TH2 depending on pt and eta (or the
674 corresponding value in case of configuration)
675 Note: In case values are outside of bin ranges, the closest bin value is used
676*/
677//______________________________________________________________________________
679 double& dEfficiencyScaleFactor, double dVars[])
680{
681 double dPt = dVars[0];
682 double dEta = dVars[1];
683
684 const TH2* hHist = dynamic_cast<const TH2*>(oObject);
685
686 if (!hHist)
687 {
688 // ATH_MSG_ERROR("Problem with casting TObject of type "<<oObject->ClassName()<<" to TH2F");
690 }
691
692 // protect values from underflow bins
693 dPt = std::max(dPt,hHist->GetXaxis()->GetXmin());
694 dEta = std::max(dEta,hHist->GetYaxis()->GetXmin());
695 // protect values from overflow bins (times .999 to keep it inside last bin)
696 dPt = std::min(dPt,hHist->GetXaxis()->GetXmax() * .999);
697 dEta = std::min(dEta,hHist->GetYaxis()->GetXmax() * .999);
698
699 // get bin from TH2 depending on x and y values; finally set the scale factor
700 int iBin = hHist->FindFixBin(dPt,dEta);
701 dEfficiencyScaleFactor = hHist->GetBinContent(iBin);
703}
704
705/*
706 Find the particular value in TF1 depending on pt and eta (or the corresponding
707 value in case of configuration)
708*/
709//______________________________________________________________________________
711 double& dEfficiencyScaleFactor, double dVars[])
712{
713 double dPt = dVars[0];
714 double dEta = dVars[1];
715
716 const TF1* fFunc = static_cast<const TF1*>(oObject);
717
718 if (!fFunc)
719 {
720 // ATH_MSG_ERROR("Problem with casting TObject of type "<<oObject->ClassName()<<" to TF1");
722 }
723
724 // evaluate TFunction and set scale factor
725 dEfficiencyScaleFactor = fFunc->Eval(dPt, dEta);
727}
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
std::string PathResolverFindCalibFile(const std::string &logical_file_name)
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.
Class to wrap a set of SystematicVariations.
std::string name() const
returns: the systematics joined into a single string.
void insert(const SystematicVariation &systematic)
description: insert a systematic into the set
size_t size() const
returns: size of the set
Helper class to provide type-safe access to aux data.
Definition Decorator.h:59
bool isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
std::string ConvertDecayModeToString(const int iDecayMode) const
std::function< double(const xAOD::TauJet &xTau)> m_fY
virtual StatusCode applySystematicVariation(const CP::SystematicSet &sSystematicSet)
configure this tool for the given list of systematic variations.
virtual CP::CorrectionCode getValue(const std::string &sHistName, const xAOD::TauJet &xTau, double &dEfficiencyScaleFactor) const
CommonEfficiencyTool(const std::string &sName)
Create a proper constructor for Athena.
std::map< std::string, tTupleObjectFunc > tSFMAP
Gaudi::Property< std::string > m_sInputFilePath
virtual CP::SystematicSet affectingSystematics() const
returns: the list of all systematics this tool can be affected by
virtual CP::SystematicSet recommendedSystematics() const
returns: the list of all systematics this tool recommends to use
static CP::CorrectionCode getValueTH1(const TObject *oObject, double &dEfficiencyScaleFactor, double dVars[])
std::tuple< TObject *, CP::CorrectionCode(*)(const TObject *oObject, double &dEfficiencyScaleFactor, double dVars[]) > tTupleObjectFunc
std::unordered_map< CP::SystematicSet, std::string > m_mSystematicSets
void addHistogramToSFMap(TKey *kKey, const std::string &sKeyName)
virtual CP::CorrectionCode getEfficiencyScaleFactor(const xAOD::TauJet &tau, double &dEfficiencyScaleFactor, unsigned int iRunNumber=0)
Declare the interface that the class provides.
static CP::CorrectionCode getValueTH2(const TObject *oObject, double &dEfficiencyScaleFactor, double dVars[])
virtual StatusCode initialize()
Dummy implementation of the initialisation function.
std::map< std::string, std::string > m_mSystematicsHistNames
virtual bool isAffectedBySystematic(const CP::SystematicVariation &systematic) const
returns: whether this tool is affected by the given systematics
std::function< double(const xAOD::TauJet &xTau)> m_fX
Gaudi::Property< std::string > m_sVarName
virtual CP::CorrectionCode applyEfficiencyScaleFactor(const xAOD::TauJet &xTau, unsigned int iRunNumber=0)
Decorate the tau with its efficiency.
static CP::CorrectionCode getValueTF1(const TObject *oObject, double &dEfficiencyScaleFactor, double dVars[])
std::string ConvertProngToString(const int iProngness) const
AsgTool(const std::string &name)
Constructor specifying the tool instance's name.
Definition AsgTool.cxx:58
bool panTauDetail(TauJetParameters::PanTauDetails panTauDetail, int &value) const
Get and set values of pantau details variables via enum.
size_t nTracks(TauJetParameters::TauTrackFlag flag=TauJetParameters::TauTrackFlag::classifiedCharged) const
double finalTauEta(const xAOD::TauJet &xTau)
return MVA based tau eta
TruthMatchedParticleType getTruthParticleType(const xAOD::TauJet &xTau)
return TauJet match type
void split(const std::string &sInput, const char cDelim, std::vector< std::string > &vOut)
double finalTauPt(const xAOD::TauJet &xTau)
return MVA based tau pt in GeV
double finalTauAbsEta(const xAOD::TauJet &xTau)
return MVA based absolute tau eta
double truthVisTauPt(const xAOD::TauJet &xTau)
return truth match visible tau pt in GeV (if hadronic truth tau match)
double truthTauAbsEta(const xAOD::TauJet &xTau)
return truth match tau eta (if hadronic truth tau match)
double truthDecayMode(const xAOD::TauJet &xTau)
return truth decay mode (if hadronic truth tau match)
TauJet_v3 TauJet
Definition of the current "tau version".