test
CMS 3D CMS Logo

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