ATLAS Offline Software
Loading...
Searching...
No Matches
semilCorr.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6
7using namespace std;
8
9semilCorr::semilCorr(const TString& fIn, const string& /*suffix*/, bool DebugIn){
10 m_Debug = DebugIn;
11 m_f.reset (TFile::Open(fIn));
12 m_etas.push_back(0);
13 m_etas.push_back(0.8);
14 m_etas.push_back(1.2);
15 m_etas.push_back(1.7);
16 m_etas.push_back(2.1);
17 m_etas.push_back(2.5);
18
19 vector<string> etastr;
20 etastr.push_back("_e0");
21 etastr.push_back("_e1");
22 etastr.push_back("_e2");
23 etastr.push_back("_e3");
24 etastr.push_back("_e4");
25
26 vector<string> prefix;
27 prefix.push_back("corr");
28 prefix.push_back("tagSyst");
29 prefix.push_back("fragSyst");
30 prefix.push_back("decaySyst");
31 prefix.push_back("msSyst");
32 prefix.push_back("idSyst");
33 prefix.push_back("decayRewSyst");
34 prefix.push_back("muRewSyst");
35 prefix.push_back("corrIncl");
36
37 for(unsigned int j = 0; j<prefix.size(); j++){
38 vector<std::unique_ptr<TH1F> > corr;
39 for(unsigned int i = 0; i<etastr.size(); i++){
40 corr.emplace_back((TH1F*) m_f->Get((prefix[j]+etastr[i]).c_str()));
41 }
42 m_histos.push_back(std::move(corr));
43 }
44}
45
46
49
50float semilCorr::getSemilCorrToIncl(const TLorentzVector& jet,
51 const TLorentzVector& mu)
52{
53 return getSemilCorrToIncl(jet,mu,m_histos[0]);
54}
55
56float semilCorr::getBjetCorrToIncl(const TLorentzVector& jet,
57 const TLorentzVector& mu)
58{
59 return getSemilCorrToIncl(jet,mu,m_histos[8]);
60}
61
62float semilCorr::getSemilCorrToIncl(const TLorentzVector& jet,
63 const TLorentzVector& mu,
64 const vector<std::unique_ptr<TH1F> >& histos)
65{
66 TLorentzVector jetmu = jet+mu;
67 //correction to get things to 1 (or to pttruth), not to reference
68 float corr = getResponse(jetmu.Pt(), jetmu.Eta(), histos);
69 return corr;
70}
71
73{
74 vector<int> indices;
75 if(syst == semilCorr::ALL){
76 indices.push_back(((int)semilCorr::TAGGINGWEIGHT));
77 indices.push_back(((int)semilCorr::FRAGMENTATION));
78 indices.push_back(((int)semilCorr::DECAY));
79 indices.push_back(((int)semilCorr::MSRESO));
80 indices.push_back(((int)semilCorr::IDRESO));
81 indices.push_back(((int)semilCorr::DECAYREW));
82 indices.push_back(((int)semilCorr::MUONSPECTRUM));
83 }else
84 indices.push_back((int) syst);
85 return indices;
86}
87
88float semilCorr::getSemilCorrToInclSyst(const TLorentzVector& jet,
89 const TLorentzVector& mu,
91{
92 //vector<int> indices = getHistoIndices(up,syst);
93 vector<int> indices1 = getHistoIndices(syst);
94 float systr = 0;
95 for(unsigned int i = 0; i<indices1.size(); i++){
96 systr += pow(getSemilCorrToIncl(jet,mu,m_histos[indices1[i]]),2);
97 }
98 return sqrt(systr);
99}
100
101float semilCorr::getResponse(float pt, float eta, const vector<std::unique_ptr<TH1F> >& h)
102{
103 float usePt = pt;
104 int histbin = -1;
105 for(unsigned int i = 0; i<m_etas.size()-1; i++){
106 if(fabs(eta)>=m_etas[i] && fabs(eta)<m_etas[i+1])
107 histbin = i;
108 }
109
110 if(histbin == -1) histbin = h.size()-1;
111 float min = h[histbin]->GetBinCenter(1);
112 float max = h[histbin]->GetBinCenter(h[histbin]->GetNbinsX());
113 if(pt>max)
114 usePt = max-0.001;
115 if(pt<min) usePt = min+0.001;
116 float result = h[histbin]->Interpolate(usePt);
117 return result;
118}
119
Scalar eta() const
pseudorapidity method
constexpr int pow(int base, int exp) noexcept
#define min(a, b)
Definition cfImp.cxx:40
#define max(a, b)
Definition cfImp.cxx:41
float getResponse(float pt, float eta, const std::vector< std::unique_ptr< TH1F > > &h)
bool m_Debug
Definition semilCorr.h:36
semilCorr(const TString &fIn, const std::string &suffix="", bool DebugIn=false)
Definition semilCorr.cxx:9
std::unique_ptr< TFile > m_f
Definition semilCorr.h:34
std::vector< float > m_etas
Definition semilCorr.h:32
float getSemilCorrToInclSyst(const TLorentzVector &jet, const TLorentzVector &mu, semilCorr::Systematics syst=semilCorr::ALL)
Definition semilCorr.cxx:88
std::vector< int > getHistoIndices(semilCorr::Systematics syst)
Definition semilCorr.cxx:72
@ MUONSPECTRUM
Definition semilCorr.h:26
@ TAGGINGWEIGHT
Definition semilCorr.h:20
@ FRAGMENTATION
Definition semilCorr.h:21
float getSemilCorrToIncl(const TLorentzVector &jet, const TLorentzVector &mu, const std::vector< std::unique_ptr< TH1F > > &histos)
Definition semilCorr.cxx:62
float getBjetCorrToIncl(const TLorentzVector &jet, const TLorentzVector &mu)
Definition semilCorr.cxx:56
std::vector< std::vector< std::unique_ptr< TH1F > > > m_histos
Definition semilCorr.h:31
STL namespace.