ATLAS Offline Software
Loading...
Searching...
No Matches
METSignificance.h
Go to the documentation of this file.
1
2/*
3 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
4*/
5// METSignificance.h
6// Header file for class METSignificance
7// Author: P. Francavilla<francav@cern.ch>
8// Author: D. Schaefer <schae@cern.ch>
10#ifndef METUTILITIES_MET_METSIGNIFICANCE_H
11#define METUTILITIES_MET_METSIGNIFICANCE_H 1
12
13// STL includes
14#include <string>
15#include <tuple>
16
17// ROOT includes
18#include "TH2F.h"
19#include "TFile.h"
20
21// FrameWork includes
22#include "AsgTools/ToolHandle.h"
23#include "AsgTools/AsgTool.h"
25
26// METInterface includes
28
29// EDM includes
31
32// Tool interfaces
38
39namespace met {
40
41 // typedefs
42 typedef ElementLink<xAOD::IParticleContainer> obj_link_t;
43
45 : public asg::AsgTool,
46 virtual public IMETSignificance
47
48 {
49 // This macro defines the constructor with the interface declaration
51
52 public:
53
54 // Copy constructor:
55
57 METSignificance(const std::string& name);
58
61
62 // Athena algtool's Hooks
63 StatusCode initialize();
64 StatusCode finalize();
65
66 StatusCode varianceMET(xAOD::MissingETContainer* metCont, float avgmu, const std::string& jetTermName, const std::string& softTermName, const std::string& totalMETName);
67
68 // Rotates the phi direction of the object resolutions & recomputes the MET significance
69 StatusCode RotateToPhi(float phi);
70
71 // Subtracts the vector lambda from the MET & recomputes the MET signficance in new MET - lambda direction
72 StatusCode SetLambda(const float px, const float py, const bool GeV=true);
73
74 double GetMETOverSqrtSumET() const { if(m_sumet>0.0) return (m_met/std::sqrt(m_sumet)); return -1.0; }
75 double GetMETOverSqrtHT () const { if(m_ht>0.0) return (m_met/std::sqrt(m_ht)); return -1.0; }
76 double GetSignificance() const { if(m_significance>0.0) return std::sqrt(m_significance); return -1.0; }
77 double GetSigDirectional() const { if(m_VarL>0.0) return m_met/std::sqrt(m_VarL); return -1.0; }
78 double GetRho() const { return m_rho; }
79 double GetVarL() const { return m_VarL; }
80 double GetVarT() const { return m_VarT; }
81 double GetTermVarL(const int term) const { if(m_term_VarL.find(term)!=m_term_VarL.end()) return m_term_VarL.find(term)->second; return -1.0e3; }
82 double GetTermVarT(const int term) const { if(m_term_VarT.find(term)!=m_term_VarT.end()) return m_term_VarT.find(term)->second; return -1.0e3; }
83
84 private:
85
88
89 ToolHandle<IJetCalibrationTool> m_jetCalibTool {this, "jetCalibTool", "", "the IJetCalibrationTool we use"};
90 ToolHandle<CP::IMuonCalibrationAndSmearingTool> m_muonCalibrationAndSmearingTool {this, "MuonCalibTool", "", "the IMuonCalibrationAndSmearingTool we use"};
91 ToolHandle<CP::IEgammaCalibrationAndSmearingTool> m_egammaCalibTool {this, "egammaCalibTool", "", "the IEgammaCalibrationAndSmearingTool we use"};
92 ToolHandle<ITauToolBase> m_tauCombinedTES {this, "tauCombinedTES", "", "the TauCombinedTES tool we use"};
93
94 Gaudi::Property<std::string> m_jetCalibConfig{this, "JetCalibConfig", "", "Config file for jet calibration/resolution"};
95 Gaudi::Property<std::string> m_jetCalibSeq{this, "JetCalibSequence", "", "Calibration sequence used for input jets"};
96 Gaudi::Property<std::string> m_jetCalibArea{this, "JetCalibArea", "", "CalibArea for input jet calibration"};
97 Gaudi::Property<int> m_muonCalibMode{this, "MuonCalibMode", -1, "CalibMode for muon momentum calibration"};
98 Gaudi::Property<std::string> m_egESModel{this, "EgammaESModel", "", "ESModel for egamma calibration"};
99 Gaudi::Property<std::string> m_egDecorrModel{this, "EgammaDecorrelationModel", "", "Decorrelation model for egamma calibration"};
100 Gaudi::Property<bool> m_egUseFastsim{this, "EgammaUseFastsim", false, "Use fastsim resolutions for egamma?"};
101 Gaudi::Property<std::string> m_tauTESConfig{this, "TauTESConfig", "CombinedTES_R22_Round2.5_v2.root", "Config file for tau energy scale calibration"};
102 Gaudi::Property<bool> m_tauUseMVARes{this, "TauUseMVAResolution", true, "Use MVA resolution for taus?"};
103
104
105 StatusCode AddMuon (const xAOD::IParticle* obj, float &pt_reso, float &phi_reso, float avgmu);
106 StatusCode AddElectron(const xAOD::IParticle* obj, float &pt_reso, float &phi_reso, float avgmu);
107 StatusCode AddPhoton (const xAOD::IParticle* obj, float &pt_reso, float &phi_reso);
108 StatusCode AddJet (const xAOD::IParticle* obj, float &pt_reso, float &phi_reso, float &avgmu);
109 StatusCode AddTau (const xAOD::IParticle* obj, float &pt_reso, float &phi_reso);
110 void AddSoftTerm(const xAOD::MissingET* soft, const TVector3 &met_vect, double (&particle_sum)[2][2]);
111
112 double GetPUProb(double jet_eta, double jet_phi,double jet_pt, double jet_jvt, double jet_fjvt, float avgmu);
113 double GetPhiUnc(double jet_eta, double jet_phi,double jet_pt);
114 static unsigned int getEtaBin(double jet_eta);
115
116 std::tuple<double,double,double> CovMatrixRotation(double var_x, double var_y, double cv_xy, double Phi);
117
118 double Significance_LT(double Numerator, double var_parall, double var_perpen, double cov);
119
120 void InvertMatrix(double (&mat)[2][2], double (&m)[2][2]);
121 void AddMatrix(double (&X)[2][2],double (&Y)[2][2], double (&mat_new)[2][2]);
122 void RotateXY(const double (&mat)[2][2], double (&mat_new)[2][2], double phi);
123
124 // soft term bias
125 double BiasPtSoftdir(const double PtSoft);
126 double VarparPtSoftdir(const double PtSoft, const double SoftSumet);
127
128 // pthard - parameterization
129 double Bias_PtSoftParall(const double PtSoft_Parall);
130 double Var_Ptsoft(const double PtSoft);
131
132 // Fill Reso map
133 void AddResoMap(const double varL, const double varT, const double CvTV, const int term);
134
135 double m_GeV;
136
141
142 TVector3 m_met_vect;
143 TVector3 m_soft_vect;
145 TVector3 m_lamda_vect;
146
153
156
157 // set limits
160
163 double m_rho;
164 double m_VarL;
165 double m_VarT;
166 double m_CvLT;
167
171
172 std::map<int, double> m_term_VarL;
173 std::map<int, double> m_term_VarT;
174 std::map<int, double> m_term_CvLT;
175
176 double m_met;
177 double m_metx;
178 double m_mety;
179 double m_metphi;
180 double m_metsoft;
182 double m_ht;
183 double m_sumet;
184
185 // Jet Uncertainties
186 TFile *m_file;
190
191 std::string m_configPrefix;
193 std::string m_JetResoAux;
194 std::string m_EMuResoAux;
195 std::string m_JetCollection;
196
197 };
198
199} //> end namespace met
200#endif //> !METUTILITIES_MET_METSIGNIFICANCE_H
#define ASG_TOOL_CLASS(CLASSNAME, INT1)
Base class for the dual-use tool implementation classes.
Definition AsgTool.h:47
virtual ~METSignificance()
Destructor:
double GetRho() const
double GetTermVarL(const int term) const
METSignificance(const std::string &name)
Constructor with parameters:
void AddSoftTerm(const xAOD::MissingET *soft, const TVector3 &met_vect, double(&particle_sum)[2][2])
double GetTermVarT(const int term) const
std::tuple< double, double, double > CovMatrixRotation(double var_x, double var_y, double cv_xy, double Phi)
double Var_Ptsoft(const double PtSoft)
ToolHandle< ITauToolBase > m_tauCombinedTES
double Significance_LT(double Numerator, double var_parall, double var_perpen, double cov)
StatusCode RotateToPhi(float phi)
double GetMETOverSqrtSumET() const
Gaudi::Property< bool > m_tauUseMVARes
Gaudi::Property< std::string > m_tauTESConfig
Gaudi::Property< std::string > m_egESModel
double BiasPtSoftdir(const double PtSoft)
Parameterization with PtSoft Direction //.
void AddMatrix(double(&X)[2][2], double(&Y)[2][2], double(&mat_new)[2][2])
Gaudi::Property< bool > m_egUseFastsim
StatusCode initialize()
Dummy implementation of the initialisation function.
Gaudi::Property< std::string > m_jetCalibSeq
std::map< int, double > m_term_VarT
ToolHandle< IJetCalibrationTool > m_jetCalibTool
double GetSigDirectional() const
double GetSignificance() const
double GetVarL() const
static unsigned int getEtaBin(double jet_eta)
std::map< int, double > m_term_VarL
void AddResoMap(const double varL, const double varT, const double CvTV, const int term)
Gaudi::Property< std::string > m_egDecorrModel
Gaudi::Property< int > m_muonCalibMode
void RotateXY(const double(&mat)[2][2], double(&mat_new)[2][2], double phi)
StatusCode AddPhoton(const xAOD::IParticle *obj, float &pt_reso, float &phi_reso)
std::string m_configJetPhiResoFile
double GetPhiUnc(double jet_eta, double jet_phi, double jet_pt)
METSignificance()
Default constructor:
StatusCode SetLambda(const float px, const float py, const bool GeV=true)
StatusCode AddElectron(const xAOD::IParticle *obj, float &pt_reso, float &phi_reso, float avgmu)
double Bias_PtSoftParall(const double PtSoft_Parall)
Gaudi::Property< std::string > m_jetCalibArea
ToolHandle< CP::IEgammaCalibrationAndSmearingTool > m_egammaCalibTool
StatusCode AddTau(const xAOD::IParticle *obj, float &pt_reso, float &phi_reso)
StatusCode varianceMET(xAOD::MissingETContainer *metCont, float avgmu, const std::string &jetTermName, const std::string &softTermName, const std::string &totalMETName)
std::map< int, double > m_term_CvLT
Gaudi::Property< std::string > m_jetCalibConfig
void InvertMatrix(double(&mat)[2][2], double(&m)[2][2])
double VarparPtSoftdir(const double PtSoft, const double SoftSumet)
double GetVarT() const
double GetMETOverSqrtHT() const
StatusCode AddMuon(const xAOD::IParticle *obj, float &pt_reso, float &phi_reso, float avgmu)
ToolHandle< CP::IMuonCalibrationAndSmearingTool > m_muonCalibrationAndSmearingTool
double GetPUProb(double jet_eta, double jet_phi, double jet_pt, double jet_jvt, double jet_fjvt, float avgmu)
StatusCode AddJet(const xAOD::IParticle *obj, float &pt_reso, float &phi_reso, float &avgmu)
Class providing the definition of the 4-vector interface.
ElementLink< xAOD::IParticleContainer > obj_link_t
MissingET_v1 MissingET
Version control by type defintion.