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 {
23  std::size_t operator()(const reco::CandidatePtr& k) const {
24  if (k.refCore().isTransient())
25  return (unsigned long)k.refCore().productPtr() ^ k.key();
26  else
27  return k.refCore().id().processIndex() ^ k.refCore().id().productIndex() ^ k.key();
28  }
29  };
30 } // namespace
31 
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 
45 
47  const std::vector<edm::Handle<reco::CandidateView> >& leptons,
48  const edm::Handle<edm::View<reco::Candidate> >& pfCandidatesH,
49  double rho,
50  JME::JetResolution& resPtObj,
51  JME::JetResolution& resPhiObj,
53  bool isRealData,
54  double& sumPtUnclustered) {
55  //pfcandidates
56  const edm::View<reco::Candidate>* pfCandidates = pfCandidatesH.product();
57 
58  // metsig covariance
59  double cov_xx = 0;
60  double cov_xy = 0;
61  double cov_yy = 0;
62 
63  // for lepton and jet subtraction
64  std::unordered_set<reco::CandidatePtr, ptr_hash> footprint;
65 
66  // subtract leptons out of sumPtUnclustered
67  for (const auto& lep_i : leptons) {
68  for (const auto& lep : *lep_i) {
69  if (lep.pt() > 10) {
70  for (unsigned int n = 0; n < lep.numberOfSourceCandidatePtrs(); n++) {
71  if (lep.sourceCandidatePtr(n).isNonnull() and lep.sourceCandidatePtr(n).isAvailable()) {
72  footprint.insert(lep.sourceCandidatePtr(n));
73  }
74  }
75  }
76  }
77  }
78  // subtract jets out of sumPtUnclustered
79  for (const auto& jet : jets) {
80  // disambiguate jets and leptons
81  if (!cleanJet(jet, leptons))
82  continue;
83  for (unsigned int n = 0; n < jet.numberOfSourceCandidatePtrs(); n++) {
84  if (jet.sourceCandidatePtr(n).isNonnull() and jet.sourceCandidatePtr(n).isAvailable()) {
85  footprint.insert(jet.sourceCandidatePtr(n));
86  }
87  }
88  }
89 
90  // calculate sumPtUnclustered
91  for (size_t i = 0; i < pfCandidates->size(); ++i) {
92  // check if candidate exists in a lepton or jet
93  bool cleancand = true;
94  if (footprint.find(pfCandidates->ptrAt(i)) == footprint.end()) {
95  //dP4 recovery
96  for (const auto& it : footprint) {
97  if ((it->p4() - (*pfCandidates)[i].p4()).Et2() < 0.000025) {
98  cleancand = false;
99  break;
100  }
101  }
102  // if not, add to sumPtUnclustered
103  if (cleancand) {
104  sumPtUnclustered += (*pfCandidates)[i].pt();
105  }
106  }
107  }
108 
109  // add jets to metsig covariance matrix and subtract them from sumPtUnclustered
110  for (const auto& jet : jets) {
111  // disambiguate jets and leptons
112  if (!cleanJet(jet, leptons))
113  continue;
114 
115  double jpt = jet.pt();
116  double jeta = jet.eta();
117  double feta = std::abs(jeta);
118  double c = jet.px() / jet.pt();
119  double s = jet.py() / jet.pt();
120 
122  parameters.setJetPt(jpt).setJetEta(jeta).setRho(rho);
123 
124  // jet energy resolutions
125  double sigmapt = resPtObj.getResolution(parameters);
126  double sigmaphi = resPhiObj.getResolution(parameters);
127  // SF not needed since is already embedded in the sigma in the dataGlobalTag
128  // double sigmaSF = isRealData ? resSFObj.getScaleFactor(parameters) : 1.0;
129 
130  // split into high-pt and low-pt sector
131  if (jpt > jetThreshold_) {
132  // high-pt jets enter into the covariance matrix via JER
133 
134  double scale = 0;
135  if (feta < jetEtas_[0])
136  scale = jetParams_[0];
137  else if (feta < jetEtas_[1])
138  scale = jetParams_[1];
139  else if (feta < jetEtas_[2])
140  scale = jetParams_[2];
141  else if (feta < jetEtas_[3])
142  scale = jetParams_[3];
143  else
144  scale = jetParams_[4];
145 
146  // double dpt = sigmaSF*scale*jpt*sigmapt;
147  double dpt = scale * jpt * sigmapt;
148  double dph = jpt * sigmaphi;
149 
150  cov_xx += dpt * dpt * c * c + dph * dph * s * s;
151  cov_xy += (dpt * dpt - dph * dph) * c * s;
152  cov_yy += dph * dph * c * c + dpt * dpt * s * s;
153 
154  } else {
155  // add the (corrected) jet to the sumPtUnclustered
156  sumPtUnclustered += jpt;
157  }
158  }
159 
160  //protection against unphysical events
161  if (sumPtUnclustered < 0)
162  sumPtUnclustered = 0;
163 
164  // add pseudo-jet to metsig covariance matrix
165  cov_xx += pjetParams_[0] * pjetParams_[0] + pjetParams_[1] * pjetParams_[1] * sumPtUnclustered;
166  cov_yy += pjetParams_[0] * pjetParams_[0] + pjetParams_[1] * pjetParams_[1] * sumPtUnclustered;
167 
168  reco::METCovMatrix cov;
169  cov(0, 0) = cov_xx;
170  cov(1, 0) = cov_xy;
171  cov(0, 1) = cov_xy;
172  cov(1, 1) = cov_yy;
173 
174  return cov;
175 }
176 
178  // covariance matrix determinant
179  double det = cov(0, 0) * cov(1, 1) - cov(0, 1) * cov(1, 0);
180 
181  // invert matrix
182  double ncov_xx = cov(1, 1) / det;
183  double ncov_xy = -cov(0, 1) / det;
184  double ncov_yy = cov(0, 0) / det;
185 
186  // product of met and inverse of covariance
187  double sig = met.px() * met.px() * ncov_xx + 2 * met.px() * met.py() * ncov_xy + met.py() * met.py() * ncov_yy;
188 
189  return sig;
190 }
191 
193  const std::vector<edm::Handle<reco::CandidateView> >& leptons) {
194  for (const auto& lep_i : leptons) {
195  for (const auto& lep : *lep_i) {
196  if (reco::deltaR2(lep, jet) < dR2match_) {
197  return false;
198  }
199  }
200  }
201  return true;
202 }
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:163
Ptr< value_type > ptrAt(size_type i) const
Base class for all types of Jets.
Definition: Jet.h:20
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:39
double px() const final
x coordinate of momentum vector
ProcessIndex productIndex() const
Definition: ProductID.h:34
size_type size() const
Definition: MET.h:41
void const * productPtr() const
Definition: RefCore.h:51
RefCore const & refCore() const
Definition: Ptr.h:167
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)
METSignificance(const edm::ParameterSet &iConfig)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
double py() const final
y coordinate of momentum vector
ProductID id() const
Definition: RefCore.h:48
bool isTransient() const
Definition: RefCore.h:103
JetParameters & setJetPt(float pt)
ProcessIndex processIndex() const
Definition: ProductID.h:33