ATLAS Offline Software
Loading...
Searching...
No Matches
CommonDiTauSmearingTool.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)
10
11// ROOT include(s)
12#include "TROOT.h"
13#include "TClass.h"
14#include "TKey.h"
15#include "TH3.h"
16
17
18using namespace TauAnalysisTools;
19
20//______________________________________________________________________________
27
28/*
29 - Find the root files with smearing inputs on eos/cvmfs using PathResolver
30 - Call further functions to process and define NP strings and so on
31 - Configure to provide nominal smearings by default
32*/
34{
35 ATH_MSG_INFO( "Initializing CommonDiTauSmearingTool" );
36
37 // only read in histograms once
38 if (m_mDTSF.empty())
39 {
40 std::string sInputFilePath = PathResolverFindCalibFile(m_sInputFilePath);
41 std::unique_ptr< TFile > fSF( TFile::Open(sInputFilePath.c_str(), "READ") );
42 if(fSF == nullptr)
43 {
44 ATH_MSG_FATAL("Could not open file " << sInputFilePath.c_str());
45 return StatusCode::FAILURE;
46 }
47 ReadInputs(fSF.get(), m_mDTSF);
48 fSF->Close();
49 }
50
52
53 // load empty systematic variation by default
54 if (applySystematicVariation(CP::SystematicSet()) != StatusCode::SUCCESS )
55 return StatusCode::FAILURE;
56
57 return StatusCode::SUCCESS;
58}
59
60/*
61 Retrieve the smearing value and if requested the values for the NP's and add
62 this stuff in quadrature. Finally apply the correction to the tau pt of the
63 non-const tau.
64*/
65//______________________________________________________________________________
67{
68 // step out here if we run on data
69 if (m_bIsData)
71
72 // check which true state is requested
74 {
76 }
77
78 double dCorrection = 1.;
79 // get standard scale factor
80 CP::CorrectionCode tmpCorrectionCode = getValue("sf",
81 xDiTau,
82 dCorrection);
83 // return correction code if histogram is not available
84 if (tmpCorrectionCode != CP::CorrectionCode::Ok)
85 return tmpCorrectionCode;
86
87 // skip further process if systematic set is empty
88 if (m_sSystematicSet->size() > 0)
89 {
90 // get uncertainties summed in quadrature
91 double dTotalSystematic2 = 0;
92 double dDirection = 0;
93 for (auto syst : *m_sSystematicSet)
94 {
95 // check if systematic is available
96 auto it = m_mSystematicsHistNames.find(syst.basename());
97
98 // get uncertainty value
99 double dUncertaintySyst = 0;
100 tmpCorrectionCode = getValue(it->second,
101 xDiTau,
102 dUncertaintySyst);
103
104 // return correction code if histogram is not available
105 if (tmpCorrectionCode != CP::CorrectionCode::Ok)
106 return tmpCorrectionCode;
107
108 // needed for up/down decision
109 dDirection = syst.parameter();
110
111 // scale uncertainty with direction, i.e. +/- n*sigma
112 dUncertaintySyst *= dDirection;
113
114 // square uncertainty and add to total uncertainty
115 dTotalSystematic2 += dUncertaintySyst * dUncertaintySyst;
116 }
117
118 // now use dDirection to use up/down uncertainty
119 dDirection = (dDirection > 0) ? +1 : -1;
120
121 // finally apply uncertainty (eff * ( 1 +/- \sum )
122 dCorrection *= 1 + dDirection * std::sqrt(dTotalSystematic2);
123 }
124
125 // finally apply correction
126 xDiTau.setP4( xDiTau.pt() * dCorrection,
127 xDiTau.eta(), xDiTau.phi(), xDiTau.m());
129}
130
131/*
132 Create a non-const copy of the passed const xDiTau object and apply the
133 correction to the non-const copy.
134 */
135//______________________________________________________________________________
137 xAOD::DiTauJet*& xDiTauCopy )
138{
139
140 // A sanity check:
141 if( xDiTauCopy )
142 {
143 ATH_MSG_WARNING( "Non-null pointer received. "
144 "There's a possible memory leak!" );
145 }
146
147 // Create a new object:
148 xDiTauCopy = new xAOD::DiTauJet();
149 xDiTauCopy->makePrivateStore( xDiTau );
150
151 // Use the other function to modify this object:
152 return applyCorrection( *xDiTauCopy );
153}
154
155//=================================PRIVATE-PART=================================
156
157template<class T>
158void CommonDiTauSmearingTool::ReadInputs(TFile* fFile, std::map<std::string, T>& mMap)
159{
160 // initialize function pointer
161 m_fX = &TruthLeadPt;
163 m_fZ = &TruthDeltaR;
164
165 TKey *kKey;
166 TIter itNext(fFile->GetListOfKeys());
167 while ((kKey = (TKey*)itNext()))
168 {
169 TClass *cClass = gROOT->GetClass(kKey->GetClassName());
170 std::string sKeyName = kKey->GetName();
171
172 if (!cClass->InheritsFrom("TH3"))
173 continue;
174
175 T tObj = (T)kKey->ReadObj();
176 tObj->SetDirectory(0);
177 mMap[sKeyName] = tObj;
178 }
179 ATH_MSG_INFO("data loaded from " << fFile->GetName());
180}
181
182//______________________________________________________________________________
184{
185 std::vector<std::string> vInputFilePath;
186 split(m_sInputFilePath,'/',vInputFilePath);
187 std::string sInputFileName = vInputFilePath.back();
188
189 // creation of basic string for all NPs, e.g. "TAUS_TRUEHADTAU_SME_TES_"
190 std::vector<std::string> vSplitInputFilePath = {};
191 split(sInputFileName,'_',vSplitInputFilePath);
192 std::string sEfficiencyType = vSplitInputFilePath.at(0);
193 std::string sTruthType = vSplitInputFilePath.at(1);
194 std::transform(sEfficiencyType.begin(), sEfficiencyType.end(), sEfficiencyType.begin(), toupper);
195 std::transform(sTruthType.begin(), sTruthType.end(), sTruthType.begin(), toupper);
196 std::string sSystematicBaseString = "TAUS_"+sTruthType+"_SME_"+sEfficiencyType+"_";
197
198 // set truth type to check for in truth matching
199 if (sTruthType=="TRUEHADTAU") m_eCheckTruth = TauAnalysisTools::TruthHadronicTau;
200 if (sTruthType=="TRUEHADDITAU") m_eCheckTruth = TauAnalysisTools::TruthHadronicDiTau;
201
202 for (auto mSF : m_mDTSF)
203 {
204 // parse for nuisance parameter in histogram name
205 std::vector<std::string> vSplitNP = {};
206 split(mSF.first,'_',vSplitNP);
207 std::string sNP = vSplitNP.at(0);
208 std::string sNPUppercase = vSplitNP.at(0);
209
210 // skip nominal scale factors
211 if (sNP == "sf") continue;
212
213 // test if NP starts with a capital letter indicating that this should be recommended
214 bool bIsRecommended = false;
215 if (isupper(sNP.at(0)))
216 bIsRecommended = true;
217
218 // make sNP uppercase and build final NP entry name
219 std::transform(sNPUppercase.begin(), sNPUppercase.end(), sNPUppercase.begin(), toupper);
220 std::string sSystematicString = sSystematicBaseString+sNPUppercase;
221
222 // add all found systematics to the AffectingSystematics
223 m_sAffectingSystematics.insert(CP::SystematicVariation (sSystematicString, 1));
224 m_sAffectingSystematics.insert(CP::SystematicVariation (sSystematicString, -1));
225 // only add found uppercase systematics to the RecommendedSystematics
226 if (bIsRecommended)
227 {
228 m_sRecommendedSystematics.insert(CP::SystematicVariation (sSystematicString, 1));
229 m_sRecommendedSystematics.insert(CP::SystematicVariation (sSystematicString, -1));
230 }
231
232 ATH_MSG_DEBUG("connected histogram base name " << sNP << " with systematic " <<sSystematicString);
233 m_mSystematicsHistNames.insert({sSystematicString,sNP});
234 }
235}
236
237//______________________________________________________________________________
239 const xAOD::DiTauJet& xDiTau,
240 double& dEfficiencyScaleFactor) const
241{
242 TH3* hHist = m_mDTSF.at(sHistName);
243 if (!hHist)
244 {
245 ATH_MSG_ERROR("Histogram with name "<<sHistName<<" was not found in input file.");
247 }
248
249 double dX = m_fX(xDiTau);
250 double dY = m_fY(xDiTau);
251 double dZ = m_fZ(xDiTau);
252
253 // protect values from underflow bins
254 dX = std::max(dX,hHist->GetXaxis()->GetXmin());
255 dY = std::max(dY,hHist->GetYaxis()->GetXmin());
256 dZ = std::max(dZ,hHist->GetZaxis()->GetXmin());
257 // protect values from overflow bins (times .999 to keep it inside last bin)
258 dX = std::min(dX,hHist->GetXaxis()->GetXmax() * .999);
259 dY = std::min(dY,hHist->GetYaxis()->GetXmax() * .999);
260 dZ = std::min(dZ,hHist->GetZaxis()->GetXmax() * .999);
261
262 int iBin = hHist->FindFixBin(dX,dY,dZ);
263 dEfficiencyScaleFactor = hHist->GetBinContent(iBin);
264
266}
267
#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)
std::string PathResolverFindCalibFile(const std::string &logical_file_name)
Return value from object correction CP tools.
@ Error
Some error happened during the object correction.
@ Ok
The correction was done successfully.
Class to wrap a set of SystematicVariations.
void makePrivateStore()
Create a new (empty) private store for this object.
double(* m_fY)(const xAOD::DiTauJet &xDiTau)
CommonDiTauSmearingTool(const std::string &sName)
Create a proper constructor for Athena.
double(* m_fX)(const xAOD::DiTauJet &xDiTau)
virtual CP::CorrectionCode applyCorrection(xAOD::DiTauJet &xDiTau) override
Declare the interface that the class provides.
virtual CP::CorrectionCode correctedCopy(const xAOD::DiTauJet &xDiTau, xAOD::DiTauJet *&xDiTauCopy) override
Create a corrected copy from a constant ditau.
double(* m_fZ)(const xAOD::DiTauJet &xDiTau)
virtual CP::CorrectionCode getValue(const std::string &sHistName, const xAOD::DiTauJet &xDiTau, double &dCorrectionFactor) const
void ReadInputs(TFile *fFile, std::map< std::string, T > &mMap)
virtual StatusCode initialize() override
Dummy implementation of the initialisation function.
virtual StatusCode applySystematicVariation(const CP::SystematicSet &sSystematicSet)
configure this tool for the given list of systematic variations.
CommonSmearingTool(const std::string &sName)
Create a proper constructor for Athena.
Gaudi::Property< bool > m_bSkipTruthMatchCheck
Gaudi::Property< std::string > m_sInputFilePath
std::map< std::string, std::string > m_mSystematicsHistNames
const CP::SystematicSet * m_sSystematicSet
void setP4(double pt, double eta, double phi, double m)
Set methods for IParticle values.
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double m() const
The invariant mass of the particle.
virtual double pt() const
The transverse momentum ( ) of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
double TruthSubleadPt(const xAOD::DiTauJet &xDiTau)
return the truth vis pT of the subleading pT matched particle.
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 TruthDeltaR(const xAOD::DiTauJet &xDiTau)
return the dR of between the leading and subleading pT matched particle.
double TruthLeadPt(const xAOD::DiTauJet &xDiTau)
return the truth vis pT of the leading pT matched particle.
DiTauJet_v1 DiTauJet
Definition of the current version.
Definition DiTauJet.h:17