CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
metsig::METSignificance Class Reference

#include <METSignificance.h>

Public Member Functions

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, edm::ValueMap< float > const *weights=nullptr)
 
 METSignificance (const edm::ParameterSet &iConfig)
 
 ~METSignificance ()
 

Static Public Member Functions

static double getSignificance (const reco::METCovMatrix &cov, const reco::MET &met)
 

Private Member Functions

bool cleanJet (const reco::Jet &jet, const std::vector< edm::Handle< reco::CandidateView > > &leptons)
 

Private Attributes

double dR2match_
 
std::vector< double > jetEtas_
 
std::vector< double > jetParams_
 
double jetThreshold_
 
std::vector< double > pjetParams_
 

Detailed Description

Definition at line 35 of file METSignificance.h.

Constructor & Destructor Documentation

◆ METSignificance()

METSignificance::METSignificance ( const edm::ParameterSet iConfig)

Definition at line 32 of file METSignificance.cc.

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 }

References dR2match_, edm::ParameterSet::getParameter(), jetEtas_, jetParams_, jetThreshold_, and pjetParams_.

◆ ~METSignificance()

METSignificance::~METSignificance ( )

Definition at line 44 of file METSignificance.cc.

44 {}

Member Function Documentation

◆ cleanJet()

bool METSignificance::cleanJet ( const reco::Jet jet,
const std::vector< edm::Handle< reco::CandidateView > > &  leptons 
)
private

Definition at line 196 of file METSignificance.cc.

197  {
198  for (const auto& lep_i : leptons) {
199  for (const auto& lep : *lep_i) {
200  if (reco::deltaR2(lep, jet) < dR2match_) {
201  return false;
202  }
203  }
204  }
205  return true;
206 }

References reco::deltaR2(), metsig::jet, and HLT_2018_cff::leptons.

◆ getCovariance()

reco::METCovMatrix METSignificance::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,
edm::ValueMap< float > const *  weights = nullptr 
)

Definition at line 46 of file METSignificance.cc.

55  {
56  //pfcandidates
57  const edm::View<reco::Candidate>& pfCandidates = *pfCandidatesH;
58 
59  // metsig covariance
60  double cov_xx = 0;
61  double cov_xy = 0;
62  double cov_yy = 0;
63 
64  // for lepton and jet subtraction
65  std::unordered_set<reco::CandidatePtr, ptr_hash> footprint;
66 
67  // subtract leptons out of sumPtUnclustered
68  for (const auto& lep_i : leptons) {
69  for (const auto& lep : lep_i->ptrs()) {
70  if (lep->pt() > 10) {
71  for (unsigned int n = 0; n < lep->numberOfSourceCandidatePtrs(); n++)
72  footprint.insert(lep->sourceCandidatePtr(n));
73  }
74  }
75  }
76 
77  std::vector<bool> cleanedJets(jets.size(), false);
78  std::transform(jets.begin(), jets.end(), cleanedJets.begin(), [this, &leptons](auto const& jet) -> bool {
79  return cleanJet(jet, leptons);
80  });
81  // subtract jets out of sumPtUnclustered
82  auto iCleaned = cleanedJets.begin();
83  for (const auto& jet : jets) {
84  // disambiguate jets and leptons
85  if (!(*iCleaned++))
86  continue;
87  for (unsigned int n = 0; n < jet.numberOfSourceCandidatePtrs(); n++) {
88  footprint.insert(jet.sourceCandidatePtr(n));
89  }
90  }
91 
92  // calculate sumPtUnclustered
93  for (size_t i = 0; i < pfCandidates.size(); ++i) {
94  // check if candidate exists in a lepton or jet
95  bool cleancand = true;
96  if (footprint.find(pfCandidates.ptrAt(i)) == footprint.end()) {
97  float weight = (weights != nullptr) ? (*weights)[pfCandidates.ptrAt(i)] : 1.0;
98  //dP4 recovery
99  for (const auto& it : footprint) {
100  if (it.isNonnull() && it.isAvailable() && (reco::deltaR2(*it, pfCandidates[i]) < 0.00000025)) {
101  cleancand = false;
102  break;
103  }
104  }
105  // if not, add to sumPtUnclustered
106  if (cleancand) {
108  }
109  }
110  }
111 
112  // add jets to metsig covariance matrix and subtract them from sumPtUnclustered
113  iCleaned = cleanedJets.begin();
114  for (const auto& jet : jets) {
115  // disambiguate jets and leptons
116  if (!(*iCleaned++))
117  continue;
118 
119  double jpt = jet.pt();
120  double jeta = jet.eta();
121  double feta = std::abs(jeta);
122  double c = jet.px() / jet.pt();
123  double s = jet.py() / jet.pt();
124 
126  parameters.setJetPt(jpt).setJetEta(jeta).setRho(rho);
127 
128  // jet energy resolutions
129  double sigmapt = resPtObj.getResolution(parameters);
130  double sigmaphi = resPhiObj.getResolution(parameters);
131  // SF not needed since is already embedded in the sigma in the dataGlobalTag
132  // double sigmaSF = isRealData ? resSFObj.getScaleFactor(parameters) : 1.0;
133 
134  // split into high-pt and low-pt sector
135  if (jpt > jetThreshold_) {
136  // high-pt jets enter into the covariance matrix via JER
137 
138  double scale = 0;
139  if (feta < jetEtas_[0])
140  scale = jetParams_[0];
141  else if (feta < jetEtas_[1])
142  scale = jetParams_[1];
143  else if (feta < jetEtas_[2])
144  scale = jetParams_[2];
145  else if (feta < jetEtas_[3])
146  scale = jetParams_[3];
147  else
148  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  // add the (corrected) jet to the sumPtUnclustered
161  }
162  }
163 
164  //protection against unphysical events
165  if (sumPtUnclustered < 0)
166  sumPtUnclustered = 0;
167 
168  // add pseudo-jet to metsig covariance matrix
169  cov_xx += pjetParams_[0] * pjetParams_[0] + pjetParams_[1] * pjetParams_[1] * sumPtUnclustered;
170  cov_yy += pjetParams_[0] * pjetParams_[0] + pjetParams_[1] * pjetParams_[1] * sumPtUnclustered;
171 
172  reco::METCovMatrix cov;
173  cov(0, 0) = cov_xx;
174  cov(1, 0) = cov_xy;
175  cov(0, 1) = cov_xy;
176  cov(1, 1) = cov_yy;
177 
178  return cov;
179 }

References funct::abs(), HltBtagPostValidation_cff::c, reco::deltaR2(), JetMETHLTOfflineSource_cfi::feta, JME::JetResolution::getResolution(), mps_fire::i, metsig::jet, METSignificanceParams_cfi::jeta, singleTopDQM_cfi::jets, HLT_2018_cff::leptons, dqmiodumpmetadata::n, zmumugammaAnalyzer_cfi::pfCandidates, alignCSCRings::s, Scenarios_cff::scale, met_cff::sumPtUnclustered, HcalDetIdTransform::transform(), mps_merge::weight, and HLT_2018_cff::weights.

Referenced by cms::PFMETProducer::getMETCovMatrix(), pat::PATMETProducer::getMETCovMatrix(), and cms::METSignificanceProducer::produce().

◆ getSignificance()

double METSignificance::getSignificance ( const reco::METCovMatrix cov,
const reco::MET met 
)
static

Definition at line 181 of file METSignificance.cc.

181  {
182  // covariance matrix determinant
183  double det = cov(0, 0) * cov(1, 1) - cov(0, 1) * cov(1, 0);
184 
185  // invert matrix
186  double ncov_xx = cov(1, 1) / det;
187  double ncov_xy = -cov(0, 1) / det;
188  double ncov_yy = cov(0, 0) / det;
189 
190  // product of met and inverse of covariance
191  double sig = met.px() * met.px() * ncov_xx + 2 * met.px() * met.py() * ncov_xy + met.py() * met.py() * ncov_yy;
192 
193  return sig;
194 }

References BTaggingMonitor_cfi::met.

Referenced by pat::PATMETProducer::produce(), CorrectedPATMETProducer::produce(), and cms::METSignificanceProducer::produce().

Member Data Documentation

◆ dR2match_

double metsig::METSignificance::dR2match_
private

Definition at line 57 of file METSignificance.h.

Referenced by METSignificance().

◆ jetEtas_

std::vector<double> metsig::METSignificance::jetEtas_
private

Definition at line 58 of file METSignificance.h.

Referenced by METSignificance().

◆ jetParams_

std::vector<double> metsig::METSignificance::jetParams_
private

Definition at line 59 of file METSignificance.h.

Referenced by METSignificance().

◆ jetThreshold_

double metsig::METSignificance::jetThreshold_
private

Definition at line 56 of file METSignificance.h.

Referenced by METSignificance().

◆ pjetParams_

std::vector<double> metsig::METSignificance::pjetParams_
private

Definition at line 60 of file METSignificance.h.

Referenced by METSignificance().

met_cff.sumPtUnclustered
sumPtUnclustered
Definition: met_cff.py:20
metsig::METSignificance::jetEtas_
std::vector< double > jetEtas_
Definition: METSignificance.h:58
zmumugammaAnalyzer_cfi.pfCandidates
pfCandidates
Definition: zmumugammaAnalyzer_cfi.py:11
HLT_2018_cff.weights
weights
Definition: HLT_2018_cff.py:87167
mps_fire.i
i
Definition: mps_fire.py:355
dqmiodumpmetadata.n
n
Definition: dqmiodumpmetadata.py:28
mps_merge.weight
weight
Definition: mps_merge.py:88
METSignificanceParams_cfi.jeta
jeta
Definition: METSignificanceParams_cfi.py:12
JME::JetParameters
Definition: JetResolutionObject.h:90
singleTopDQM_cfi.jets
jets
Definition: singleTopDQM_cfi.py:42
BTaggingMonitor_cfi.met
met
Definition: BTaggingMonitor_cfi.py:84
parameters
parameters
Definition: BeamSpot_PayloadInspector.cc:14
metsig::METSignificance::pjetParams_
std::vector< double > pjetParams_
Definition: METSignificance.h:60
alignCSCRings.s
s
Definition: alignCSCRings.py:92
jpt
Definition: JetPlusTrackCorrector.h:37
HcalDetIdTransform::transform
unsigned transform(const HcalDetId &id, unsigned transformCode)
Definition: HcalDetIdTransform.cc:7
DDAxes::rho
edm::View
Definition: CaloClusterFwd.h:14
Scenarios_cff.scale
scale
Definition: Scenarios_cff.py:2186
JetMETHLTOfflineSource_cfi.feta
feta
Definition: JetMETHLTOfflineSource_cfi.py:30
edm::ParameterSet
Definition: ParameterSet.h:36
reco::deltaR2
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
HLT_2018_cff.leptons
leptons
Definition: HLT_2018_cff.py:24820
metsig::METSignificance::dR2match_
double dR2match_
Definition: METSignificance.h:57
metsig::METSignificance::jetThreshold_
double jetThreshold_
Definition: METSignificance.h:56
HltBtagPostValidation_cff.c
c
Definition: HltBtagPostValidation_cff.py:31
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
metsig::METSignificance::jetParams_
std::vector< double > jetParams_
Definition: METSignificance.h:59
metsig::jet
Definition: SignAlgoResolutions.h:47
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
metsig::METSignificance::cleanJet
bool cleanJet(const reco::Jet &jet, const std::vector< edm::Handle< reco::CandidateView > > &leptons)
Definition: METSignificance.cc:196
weight
Definition: weight.py:1
reco::METCovMatrix
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:39
JME::JetResolution::getResolution
float getResolution(const JetParameters &parameters) const
Definition: JetResolution.cc:31