CMS 3D CMS Logo

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