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  useDeltaRforFootprint_ = cfgParams.getParameter<bool>("useDeltaRforFootprint");
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  double& sumPtUnclustered) {
59 
60  //pfcandidates
61  const edm::View<reco::Candidate>* pfCandidates=pfCandidatesH.product();
62 
63  // metsig covariance
64  double cov_xx = 0;
65  double cov_xy = 0;
66  double cov_yy = 0;
67 
68  // for lepton and jet subtraction
69  std::unordered_set<reco::CandidatePtr,ptr_hash> footprint;
70 
71  // subtract leptons out of sumPtUnclustered
72  for ( const auto& lep_i : leptons ) {
73  for (const auto& lep : lep_i->ptrs()) {
74  if (lep->pt() > 10) {
75  for (unsigned int n = 0; n < lep->numberOfSourceCandidatePtrs(); n++)
76  footprint.insert(lep->sourceCandidatePtr(n));
77  }
78  }
79  }
80  // subtract jets out of sumPtUnclustered
81  for(const auto& jet : jets) {
82 
83  // disambiguate jets and leptons
84  if(!cleanJet(jet, leptons) ) continue;
85  for( unsigned int n=0; n < jet.numberOfSourceCandidatePtrs(); n++){
86  footprint.insert(jet.sourceCandidatePtr(n));
87  }
88 
89  }
90 
91  // calculate sumPtUnclustered
92  for(size_t i = 0; i< pfCandidates->size(); ++i) {
93 
94  // check if candidate exists in a lepton or jet
95  bool cleancand = true;
96  if(footprint.find( pfCandidates->ptrAt(i) )==footprint.end()) {
97 
98  //dP4 recovery
99  for( const auto& it : footprint) {
100  // Special treatment for PUPPI with dR, since jet candidates (puppi) and MET candidates (puppiForMet)
101  // can't be matched through the sourceCandidatePtrs and may have different energy, but same direction.
102  if((it.isNonnull()) && (it.isAvailable()) &&
103  (((!useDeltaRforFootprint_) && ((it->p4()-(*pfCandidates)[i].p4()).Et2()<0.000025)) ||
104  (( useDeltaRforFootprint_) && (reco::deltaR2(*it,(*pfCandidates)[i])<0.00000025)))){
105 
106  cleancand = false;
107  break;
108  }
109  }
110  // if not, add to sumPtUnclustered
111  if( cleancand ){
112  sumPtUnclustered += (*pfCandidates)[i].pt();
113  }
114  }
115  }
116 
117  // add jets to metsig covariance matrix and subtract them from sumPtUnclustered
118  for(const auto& jet : jets) {
119 
120  // disambiguate jets and leptons
121  if(!cleanJet(jet, leptons) ) continue;
122 
123  double jpt = jet.pt();
124  double jeta = jet.eta();
125  double feta = std::abs(jeta);
126  double c = jet.px()/jet.pt();
127  double s = jet.py()/jet.pt();
128 
130  parameters.setJetPt(jpt).setJetEta(jeta).setRho(rho);
131 
132  // jet energy resolutions
133  double sigmapt = resPtObj.getResolution(parameters);
134  double sigmaphi = resPhiObj.getResolution(parameters);
135  // SF not needed since is already embedded in the sigma in the dataGlobalTag
136  // double sigmaSF = isRealData ? resSFObj.getScaleFactor(parameters) : 1.0;
137 
138  // split into high-pt and low-pt sector
139  if( jpt > jetThreshold_ ){
140  // high-pt jets enter into the covariance matrix via JER
141 
142  double scale = 0;
143  if(feta<jetEtas_[0]) scale = jetParams_[0];
144  else if(feta<jetEtas_[1]) scale = jetParams_[1];
145  else if(feta<jetEtas_[2]) scale = jetParams_[2];
146  else if(feta<jetEtas_[3]) scale = jetParams_[3];
147  else scale = jetParams_[4];
148 
149  // double dpt = sigmaSF*scale*jpt*sigmapt;
150  double dpt = scale*jpt*sigmapt;
151  double dph = jpt*sigmaphi;
152 
153  cov_xx += dpt*dpt*c*c + dph*dph*s*s;
154  cov_xy += (dpt*dpt-dph*dph)*c*s;
155  cov_yy += dph*dph*c*c + dpt*dpt*s*s;
156 
157  } else {
158 
159  // add the (corrected) jet to the sumPtUnclustered
160  sumPtUnclustered += jpt;
161 
162  }
163 
164  }
165 
166  //protection against unphysical events
167  if(sumPtUnclustered<0) sumPtUnclustered=0;
168 
169  // add pseudo-jet to metsig covariance matrix
170  cov_xx += pjetParams_[0]*pjetParams_[0] + pjetParams_[1]*pjetParams_[1]*sumPtUnclustered;
171  cov_yy += pjetParams_[0]*pjetParams_[0] + pjetParams_[1]*pjetParams_[1]*sumPtUnclustered;
172 
173  reco::METCovMatrix cov;
174  cov(0,0) = cov_xx;
175  cov(1,0) = cov_xy;
176  cov(0,1) = cov_xy;
177  cov(1,1) = cov_yy;
178 
179  return cov;
180 }
181 
182 double
184 
185  // covariance matrix determinant
186  double det = cov(0,0)*cov(1,1) - cov(0,1)*cov(1,0);
187 
188  // invert matrix
189  double ncov_xx = cov(1,1) / det;
190  double ncov_xy = -cov(0,1) / det;
191  double ncov_yy = cov(0,0) / det;
192 
193  // product of met and inverse of covariance
194  double sig = met.px()*met.px()*ncov_xx + 2*met.px()*met.py()*ncov_xy + met.py()*met.py()*ncov_yy;
195 
196  return sig;
197 }
198 
199 bool
201  const std::vector< edm::Handle<reco::CandidateView> >& leptons ){
202 
203  for ( const auto& lep_i : leptons ) {
204  for( const auto& lep : *lep_i ) {
205  if ( reco::deltaR2(lep, jet) < dR2match_ ) {
206  return false;
207  }
208  }
209  }
210  return true;
211 }
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
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
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, double &sumPtUnclustered)
int k[5][pyjets_maxn]
METSignificance(const edm::ParameterSet &iConfig)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
met
===> hadronic RAZOR
double py() const final
y coordinate of momentum vector
ProductID id() const
Definition: RefCore.h:49
bool isTransient() const
Definition: RefCore.h:106
sumPtUnclustered
Definition: met_cff.py:20
JetParameters & setJetPt(float pt)
ProcessIndex processIndex() const
Definition: ProductID.h:36