-
Notifications
You must be signed in to change notification settings - Fork 4.2k
/
METSignificance.cc
212 lines (167 loc) · 6.12 KB
/
METSignificance.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
// -*- C++ -*-
//
// Package: METAlgorithms
// Class: METSignificance
//
/**\class METSignificance METSignificance.cc RecoMET/METAlgorithms/src/METSignificance.cc
Description: [one line class summary]
Implementation:
[Notes on implementation]
*/
//
// Original Author: Nathan Mirman (Cornell University)
// Created: Thu May 30 16:39:52 CDT 2013
//
//
#include "RecoMET/METAlgorithms/interface/METSignificance.h"
#include <unordered_set>
namespace {
struct ptr_hash : public std::unary_function<reco::CandidatePtr, std::size_t> {
std::size_t operator()(const reco::CandidatePtr& k) const
{
if(k.refCore().isTransient()) return (unsigned long)k.refCore().productPtr() ^ k.key();
else return k.refCore().id().processIndex() ^ k.refCore().id().productIndex() ^ k.key();
}
};
}
metsig::METSignificance::METSignificance(const edm::ParameterSet& iConfig) {
edm::ParameterSet cfgParams = iConfig.getParameter<edm::ParameterSet>("parameters");
double dRmatch = cfgParams.getParameter<double>("dRMatch");
dR2match_ = dRmatch*dRmatch;
jetThreshold_ = cfgParams.getParameter<double>("jetThreshold");
jetEtas_ = cfgParams.getParameter<std::vector<double> >("jeta");
jetParams_ = cfgParams.getParameter<std::vector<double> >("jpar");
pjetParams_ = cfgParams.getParameter<std::vector<double> >("pjpar");
}
metsig::METSignificance::~METSignificance() {
}
reco::METCovMatrix
metsig::METSignificance::getCovariance(const edm::View<reco::Jet>& jets,
const std::vector< edm::Handle<reco::CandidateView> >& leptons,
const edm::Handle<edm::View<reco::Candidate> >& pfCandidatesH,
double rho,
JME::JetResolution& resPtObj,
JME::JetResolution& resPhiObj,
JME::JetResolutionScaleFactor& resSFObj,
bool isRealData) {
//pfcandidates
const edm::View<reco::Candidate>* pfCandidates=pfCandidatesH.product();
// metsig covariance
double cov_xx = 0;
double cov_xy = 0;
double cov_yy = 0;
// for lepton and jet subtraction
std::unordered_set<reco::CandidatePtr,ptr_hash> footprint;
// subtract leptons out of sumPt
for ( const auto& lep_i : leptons ) {
for( const auto& lep : *lep_i ) {
if( lep.pt() > 10 ){
for( unsigned int n=0; n < lep.numberOfSourceCandidatePtrs(); n++ ){
if( lep.sourceCandidatePtr(n).isNonnull() and lep.sourceCandidatePtr(n).isAvailable() ){
footprint.insert(lep.sourceCandidatePtr(n));
}
}
}
}
}
// subtract jets out of sumPt
for(const auto& jet : jets) {
// disambiguate jets and leptons
if(!cleanJet(jet, leptons) ) continue;
for( unsigned int n=0; n < jet.numberOfSourceCandidatePtrs(); n++){
if( jet.sourceCandidatePtr(n).isNonnull() and jet.sourceCandidatePtr(n).isAvailable() ){
footprint.insert(jet.sourceCandidatePtr(n));
}
}
}
// calculate sumPt
double sumPt = 0;
for(size_t i = 0; i< pfCandidates->size(); ++i) {
// check if candidate exists in a lepton or jet
bool cleancand = true;
if(footprint.find( pfCandidates->ptrAt(i) )==footprint.end()) {
//dP4 recovery
for( const auto& it : footprint) {
if( (it->p4()-(*pfCandidates)[i].p4()).Et2()<0.000025 ){
cleancand = false;
break;
}
}
// if not, add to sumPt
if( cleancand ){
sumPt += (*pfCandidates)[i].pt();
}
}
}
// add jets to metsig covariance matrix and subtract them from sumPt
for(const auto& jet : jets) {
// disambiguate jets and leptons
if(!cleanJet(jet, leptons) ) continue;
double jpt = jet.pt();
double jeta = jet.eta();
double feta = std::abs(jeta);
double c = jet.px()/jet.pt();
double s = jet.py()/jet.pt();
JME::JetParameters parameters;
parameters.setJetPt(jpt).setJetEta(jeta).setRho(rho);
// jet energy resolutions
double sigmapt = resPtObj.getResolution(parameters);
double sigmaphi = resPhiObj.getResolution(parameters);
// SF not needed since is already embedded in the sigma in the dataGlobalTag
// double sigmaSF = isRealData ? resSFObj.getScaleFactor(parameters) : 1.0;
// split into high-pt and low-pt sector
if( jpt > jetThreshold_ ){
// high-pt jets enter into the covariance matrix via JER
double scale = 0;
if(feta<jetEtas_[0]) scale = jetParams_[0];
else if(feta<jetEtas_[1]) scale = jetParams_[1];
else if(feta<jetEtas_[2]) scale = jetParams_[2];
else if(feta<jetEtas_[3]) scale = jetParams_[3];
else scale = jetParams_[4];
// double dpt = sigmaSF*scale*jpt*sigmapt;
double dpt = scale*jpt*sigmapt;
double dph = jpt*sigmaphi;
cov_xx += dpt*dpt*c*c + dph*dph*s*s;
cov_xy += (dpt*dpt-dph*dph)*c*s;
cov_yy += dph*dph*c*c + dpt*dpt*s*s;
} else {
// add the (corrected) jet to the sumPt
sumPt += jpt;
}
}
//protection against unphysical events
if(sumPt<0) sumPt=0;
// add pseudo-jet to metsig covariance matrix
cov_xx += pjetParams_[0]*pjetParams_[0] + pjetParams_[1]*pjetParams_[1]*sumPt;
cov_yy += pjetParams_[0]*pjetParams_[0] + pjetParams_[1]*pjetParams_[1]*sumPt;
reco::METCovMatrix cov;
cov(0,0) = cov_xx;
cov(1,0) = cov_xy;
cov(0,1) = cov_xy;
cov(1,1) = cov_yy;
return cov;
}
double
metsig::METSignificance::getSignificance(const reco::METCovMatrix& cov, const reco::MET& met) {
// covariance matrix determinant
double det = cov(0,0)*cov(1,1) - cov(0,1)*cov(1,0);
// invert matrix
double ncov_xx = cov(1,1) / det;
double ncov_xy = -cov(0,1) / det;
double ncov_yy = cov(0,0) / det;
// product of met and inverse of covariance
double sig = met.px()*met.px()*ncov_xx + 2*met.px()*met.py()*ncov_xy + met.py()*met.py()*ncov_yy;
return sig;
}
bool
metsig::METSignificance::cleanJet(const reco::Jet& jet,
const std::vector< edm::Handle<reco::CandidateView> >& leptons ){
for ( const auto& lep_i : leptons ) {
for( const auto& lep : *lep_i ) {
if ( reco::deltaR2(lep, jet) < dR2match_ ) {
return false;
}
}
}
return true;
}