CMS 3D CMS Logo

PFMETProducerMVA.cc
Go to the documentation of this file.
2 
3 using namespace reco;
4 
5 const double dR2Min = 0.01*0.01;
6 const double dR2Max = 0.5*0.5;
7 const double dPtMatch = 0.1;
8 
10  : mvaMEtAlgo_(cfg),
11  mvaMEtAlgo_isInitialized_(false)
12 {
13  srcCorrJets_ = consumes<reco::PFJetCollection>(cfg.getParameter<edm::InputTag>("srcCorrJets"));
14  srcUncorrJets_ = consumes<reco::PFJetCollection>(cfg.getParameter<edm::InputTag>("srcUncorrJets"));
15  srcJetIds_ = consumes<edm::ValueMap<float> >(cfg.getParameter<edm::InputTag>("srcMVAPileupJetId"));
16  srcPFCandidatesView_ = consumes<reco::CandidateView>(cfg.getParameter<edm::InputTag>("srcPFCandidates"));
17  srcVertices_ = consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>("srcVertices"));
18  vInputTag srcLeptonsTags = cfg.getParameter<vInputTag>("srcLeptons");
19  for(vInputTag::const_iterator it=srcLeptonsTags.begin();it!=srcLeptonsTags.end();it++) {
20  srcLeptons_.push_back( consumes<reco::CandidateView >( *it ) );
21  }
22  mJetCorrector_ = consumes<reco::JetCorrector>(cfg.getParameter<edm::InputTag>("corrector"));
23  minNumLeptons_ = cfg.getParameter<int>("minNumLeptons");
24 
25  globalThreshold_ = cfg.getParameter<double>("globalThreshold");
26 
27  minCorrJetPt_ = cfg.getParameter<double> ("minCorrJetPt");
28  useType1_ = cfg.getParameter<bool> ("useType1");
29 
30  verbosity_ = ( cfg.exists("verbosity") ) ?
31  cfg.getParameter<int>("verbosity") : 0;
32 
33  produces<reco::PFMETCollection>();
34 }
35 
37 
39 {
40  // CV: check if the event is to be skipped
41  if ( minNumLeptons_ > 0 ) {
42  int numLeptons = 0;
43  for ( std::vector<edm::EDGetTokenT<reco::CandidateView> >::const_iterator srcLeptons_i = srcLeptons_.begin();
44  srcLeptons_i != srcLeptons_.end(); ++srcLeptons_i ) {
46  evt.getByToken(*srcLeptons_i, leptons);
47  numLeptons += leptons->size();
48  }
49  if ( !(numLeptons >= minNumLeptons_) ) {
50  LogDebug( "produce" )
51  << "<PFMETProducerMVA::produce>:" << std::endl
52  << "Run: " << evt.id().run() << ", LS: " << evt.luminosityBlock() << ", Event: " << evt.id().event() << std::endl
53  << " numLeptons = " << numLeptons << ", minNumLeptons = " << minNumLeptons_ << " --> skipping !!" << std::endl;
54 
55  reco::PFMET pfMEt;
56  auto pfMEtCollection = std::make_unique<reco::PFMETCollection>();
57  pfMEtCollection->push_back(pfMEt);
58  evt.put(std::move(pfMEtCollection));
59  return;
60  }
61  }
62 
63  //get jet IDs
65  evt.getByToken(srcJetIds_, jetIds);
66 
67  // get jets (corrected and uncorrected)
69  evt.getByToken(srcCorrJets_, corrJets);
70 
72  evt.getByToken(srcUncorrJets_, uncorrJets);
73 
75  if( useType1_ )
76  {
77  evt.getByToken(mJetCorrector_, corrector);
78  }
79 
80  edm::Handle<reco::CandidateView> pfCandidates_view;
81  evt.getByToken(srcPFCandidatesView_, pfCandidates_view);
82 
83  // get vertices
85  evt.getByToken(srcVertices_, vertices);
86  // take vertex with highest sum(trackPt) as the vertex of the "hard scatter" interaction
87  // (= first entry in vertex collection)
88  const reco::Vertex* hardScatterVertex = ( !vertices->empty() ) ?
89  &(vertices->front()) : nullptr;
90 
91  // get leptons
92  // (excluded from sum over PFCandidates when computing hadronic recoil)
93  int lId = 0;
94  bool lHasPhotons = false;
95  std::vector<reco::PUSubMETCandInfo> leptonInfo = computeLeptonInfo(srcLeptons_,*pfCandidates_view,hardScatterVertex, lId, lHasPhotons,
96  evt);
97 
98  // initialize MVA MET algorithm
99  // (this will load the BDTs, stored as GBRForrest objects;
100  // either in input ROOT files or in SQL-lite files/the Conditions Database)
101  if ( !mvaMEtAlgo_isInitialized_ ) {
104  }
105 
106  // reconstruct "standard" particle-flow missing Et
107  CommonMETData pfMEt_data = metAlgo_.run( (*pfCandidates_view), globalThreshold_);
108  SpecificPFMETData specificPfMET = pfMEtSpecificAlgo_.run( (*pfCandidates_view) );
109  const reco::Candidate::LorentzVector p4( pfMEt_data.mex, pfMEt_data.mey, 0.0, pfMEt_data.met);
110  const reco::Candidate::Point vtx(0.0, 0.0, 0.0 );
111  reco::PFMET pfMEt(specificPfMET,pfMEt_data.sumet, p4, vtx);
112  reco::Candidate::LorentzVector pfMEtP4_original = pfMEt.p4();
113 
114  // compute objects specific to MVA based MET reconstruction
115  std::vector<reco::PUSubMETCandInfo> pfCandidateInfo = computePFCandidateInfo(*pfCandidates_view, hardScatterVertex);
116  std::vector<reco::PUSubMETCandInfo> jetInfo = computeJetInfo(*uncorrJets, corrJets, *jetIds, *vertices, hardScatterVertex, *corrector,evt,es,leptonInfo,pfCandidateInfo);
117  std::vector<reco::Vertex::Point> vertexInfo = computeVertexInfo(*vertices);
118  // compute MVA based MET and estimate of its uncertainty
119  mvaMEtAlgo_.setInput(leptonInfo, jetInfo, pfCandidateInfo, vertexInfo);
120  mvaMEtAlgo_.setHasPhotons(lHasPhotons);
122  pfMEt.setP4(mvaMEtAlgo_.getMEt());
124 
125  LogDebug("produce")
126  << "Run: " << evt.id().run() << ", LS: " << evt.luminosityBlock() << ", Event: " << evt.id().event() << std::endl
127  << " PFMET: Pt = " << pfMEtP4_original.pt() << ", phi = " << pfMEtP4_original.phi() << " "
128  << "(Px = " << pfMEtP4_original.px() << ", Py = " << pfMEtP4_original.py() << ")" << std::endl
129  << " MVA MET: Pt = " << pfMEt.pt() << " phi = " << pfMEt.phi() << " (Px = " << pfMEt.px() << ", Py = " << pfMEt.py() << ")" << std::endl
130  << " Cov:" << std::endl
131  <<(mvaMEtAlgo_.getMEtCov())(0,0)<<" "<<(mvaMEtAlgo_.getMEtCov())(0,1)<<std::endl
132  <<(mvaMEtAlgo_.getMEtCov())(1,0)<<" "<<(mvaMEtAlgo_.getMEtCov())(1,1)<<std::endl << std::endl;
133 
134  // add PFMET object to the event
135  auto pfMEtCollection = std::make_unique<reco::PFMETCollection>();
136  pfMEtCollection->push_back(pfMEt);
137  evt.put(std::move(pfMEtCollection));
138 }
139 
140 std::vector<reco::PUSubMETCandInfo>
142  const reco::CandidateView& pfCandidates_view,
143  const reco::Vertex* hardScatterVertex,
144  int& lId, bool& lHasPhotons, edm::Event& evt ) {
145 
146  std::vector<reco::PUSubMETCandInfo> leptonInfo;
147 
148  for ( std::vector<edm::EDGetTokenT<reco::CandidateView > >::const_iterator srcLeptons_i = srcLeptons_.begin();
149  srcLeptons_i != srcLeptons_.end(); ++srcLeptons_i ) {
151  evt.getByToken(*srcLeptons_i, leptons);
152  for ( reco::CandidateView::const_iterator lepton1 = leptons->begin();
153  lepton1 != leptons->end(); ++lepton1 ) {
154  bool pMatch = false;
155  for ( std::vector<edm::EDGetTokenT<reco::CandidateView> >::const_iterator srcLeptons_j = srcLeptons_.begin();
156  srcLeptons_j != srcLeptons_.end(); ++srcLeptons_j ) {
158  evt.getByToken(*srcLeptons_j, leptons2);
159  for ( reco::CandidateView::const_iterator lepton2 = leptons2->begin();
160  lepton2 != leptons2->end(); ++lepton2 ) {
161  if(&(*lepton1) == &(*lepton2)) { continue; }
162  if(deltaR2(lepton1->p4(),lepton2->p4()) < dR2Max) { pMatch = true; }
163  if(pMatch && !istau(&(*lepton1)) && istau(&(*lepton2))) { pMatch = false; }
164  if(pMatch && ( (istau(&(*lepton1)) && istau(&(*lepton2))) || (!istau(&(*lepton1)) && !istau(&(*lepton2))))
165  && lepton1->pt() > lepton2->pt()) { pMatch = false; }
166  if(pMatch && lepton1->pt() == lepton2->pt()) {
167  pMatch = false;
168  for(unsigned int i0 = 0; i0 < leptonInfo.size(); i0++) {
169  if(std::abs(lepton1->pt() - leptonInfo[i0].p4().pt()) < dPtMatch) { pMatch = true; break; }
170  }
171  }
172  if(pMatch) break;
173  }
174  if(pMatch) break;
175  }
176  if(pMatch) continue;
177  reco::PUSubMETCandInfo pLeptonInfo;
178  pLeptonInfo.setP4( lepton1->p4() );
179  pLeptonInfo.setChargedEnFrac( chargedEnFrac(&(*lepton1),pfCandidates_view,hardScatterVertex) );
180  leptonInfo.push_back(pLeptonInfo);
181  if(lepton1->isPhoton()) { lHasPhotons = true; }
182  }
183  lId++;
184  }
185 
186  return leptonInfo;
187 }
188 
189 
190 std::vector<reco::PUSubMETCandInfo>
192  const edm::Handle<reco::PFJetCollection>& corrJets,
193  const edm::ValueMap<float>& jetIds,
195  const reco::Vertex* hardScatterVertex,
196  const reco::JetCorrector &iCorrector,edm::Event &iEvent,const edm::EventSetup &iSetup,
197  std::vector<reco::PUSubMETCandInfo> &iLeptons,std::vector<reco::PUSubMETCandInfo> &iCands)
198 {
199  std::vector<reco::PUSubMETCandInfo> retVal;
200  for ( reco::PFJetCollection::const_iterator uncorrJet = uncorrJets.begin();
201  uncorrJet != uncorrJets.end(); ++uncorrJet ) {
202  // for ( reco::PFJetCollection::const_iterator corrJet = corrJets.begin();
203  // corrJet != corrJets.end(); ++corrJet ) {
204  auto corrJet = corrJets->begin();
205  for( size_t cjIdx=0;cjIdx<corrJets->size();++cjIdx, ++corrJet) {
206  reco::PFJetRef corrJetRef( corrJets, cjIdx );
207 
208  // match corrected and uncorrected jets
209  if ( uncorrJet->jetArea() != corrJet->jetArea() ) continue;
210  if ( deltaR2(corrJet->p4(),uncorrJet->p4()) > dR2Min ) continue;
211 
212  // check that jet passes loose PFJet id.
213  if(!passPFLooseId(&(*uncorrJet))) continue;
214 
215  // compute jet energy correction factor
216  // (= ratio of corrected/uncorrected jet Pt)
217  //double jetEnCorrFactor = corrJet->pt()/uncorrJet->pt();
218  reco::PUSubMETCandInfo jetInfo;
219 
220  // PH: apply jet energy corrections for all Jets ignoring recommendations
221  jetInfo.setP4( corrJet->p4() );
222  double lType1Corr = 0;
223  if(useType1_) { //Compute the type 1 correction ===> This code is crap
224  double pCorr = iCorrector.correction(*uncorrJet);
225  lType1Corr = std::abs(corrJet->pt()-pCorr*uncorrJet->pt());
226  TLorentzVector pVec; pVec.SetPtEtaPhiM(lType1Corr,0,corrJet->phi(),0);
227  reco::Candidate::LorentzVector pType1Corr; pType1Corr.SetCoordinates(pVec.Px(),pVec.Py(),pVec.Pz(),pVec.E());
228  //Filter to leptons
229  bool pOnLepton = false;
230  for(unsigned int i0 = 0; i0 < iLeptons.size(); i0++) {
231  if(deltaR2(iLeptons[i0].p4(),corrJet->p4()) < dR2Max) {pOnLepton = true; break;}
232  }
233  //Add it to PF Collection
234  if(corrJet->pt() > 10 && !pOnLepton) {
235  reco::PUSubMETCandInfo pfCandidateInfo;
236  pfCandidateInfo.setP4( pType1Corr );
237  pfCandidateInfo.setDZ( -999 );
238  iCands.push_back(pfCandidateInfo);
239  }
240  //Scale
241  lType1Corr = (pCorr*uncorrJet->pt()-uncorrJet->pt());
242  lType1Corr /=corrJet->pt();
243  }
244 
245  // check that jet Pt used to compute MVA based jet id. is above threshold
246  if ( !(jetInfo.p4().pt() > minCorrJetPt_) ) continue;
247 
248 
249  jetInfo.setMvaVal( jetIds[ corrJetRef ] );
250  float chEnF = (uncorrJet->chargedEmEnergy() + uncorrJet->chargedHadronEnergy() + uncorrJet->chargedMuEnergy() )/uncorrJet->energy();
251  if(useType1_) chEnF += lType1Corr*(1-jetInfo.chargedEnFrac() );
252  jetInfo.setChargedEnFrac( chEnF );
253  retVal.push_back(jetInfo);
254  break;
255  }
256  }
257 
258  //order jets per pt
259  std::sort( retVal.begin(), retVal.end() );
260 
261  return retVal;
262 }
263 
265  const reco::Vertex* hardScatterVertex)
266 {
267  std::vector<reco::PUSubMETCandInfo> retVal;
268  for ( reco::CandidateView::const_iterator pfCandidate = pfCandidates.begin();
269  pfCandidate != pfCandidates.end(); ++pfCandidate ) {
270  double dZ = -999.; // PH: If no vertex is reconstructed in the event
271  // or PFCandidate has no track, set dZ to -999
272  if ( hardScatterVertex ) {
273  const reco::PFCandidate* pfc = dynamic_cast<const reco::PFCandidate* >( &(*pfCandidate) );
274  if( pfc != nullptr ) { //PF candidate for RECO and PAT levels
275  if ( pfc->trackRef().isNonnull() ) dZ = std::abs(pfc->trackRef()->dz(hardScatterVertex->position()));
276  else if ( pfc->gsfTrackRef().isNonnull() ) dZ = std::abs(pfc->gsfTrackRef()->dz(hardScatterVertex->position()));
277  }
278  else { //if not, then packedCandidate for miniAOD level
279  const pat::PackedCandidate* pfc = dynamic_cast<const pat::PackedCandidate* >( &(*pfCandidate) );
280  dZ = std::abs( pfc->dz( hardScatterVertex->position() ) );
281  //exact dz=zero corresponds to the -999 case for pfcandidate
282  if(dZ==0) {dZ=-999;}
283  }
284  }
285  reco::PUSubMETCandInfo pfCandidateInfo;
286  pfCandidateInfo.setP4( pfCandidate->p4() );
287  pfCandidateInfo.setDZ( dZ );
288  retVal.push_back(pfCandidateInfo);
289  }
290  return retVal;
291 }
292 
294 {
295  std::vector<reco::Vertex::Point> retVal;
296  for ( reco::VertexCollection::const_iterator vertex = vertices.begin();
297  vertex != vertices.end(); ++vertex ) {
298  if(std::abs(vertex->z()) > 24.) continue;
299  if(vertex->ndof() < 4.) continue;
300  if(vertex->position().Rho() > 2.) continue;
301  retVal.push_back(vertex->position());
302  }
303  return retVal;
304 }
306  const reco::CandidateView& pfCandidates,const reco::Vertex* hardScatterVertex) {
307  if(iCand->isMuon()) {
308  return 1;
309  }
310  if(iCand->isElectron()) {
311  return 1.;
312  }
313  if(iCand->isPhoton() ) {return chargedFracInCone(iCand, pfCandidates,hardScatterVertex);}
314  double lPtTot = 0; double lPtCharged = 0;
315  const reco::PFTau *lPFTau = nullptr;
316  lPFTau = dynamic_cast<const reco::PFTau*>(iCand);
317  if(lPFTau != nullptr) {
318  for (UInt_t i0 = 0; i0 < lPFTau->signalPFCands().size(); i0++) {
319  lPtTot += (lPFTau->signalPFCands())[i0]->pt();
320  if((lPFTau->signalPFCands())[i0]->charge() == 0) continue;
321  lPtCharged += (lPFTau->signalPFCands())[i0]->pt();
322  }
323  }
324  else {
325  const pat::Tau *lPatPFTau = nullptr;
326  lPatPFTau = dynamic_cast<const pat::Tau*>(iCand);
327  if(lPatPFTau != nullptr) {
328  for (UInt_t i0 = 0; i0 < lPatPFTau->signalCands().size(); i0++) {
329  lPtTot += (lPatPFTau->signalCands())[i0]->pt();
330  if((lPatPFTau->signalCands())[i0]->charge() == 0) continue;
331  lPtCharged += (lPatPFTau->signalCands())[i0]->pt();
332  }
333  }
334  }
335  if(lPtTot == 0) lPtTot = 1.;
336  return lPtCharged/lPtTot;
337 }
338 //Return tau id by process of elimination
340  if(iCand->isMuon()) return false;
341  if(iCand->isElectron()) return false;
342  if(iCand->isPhoton()) return false;
343  return true;
344 }
346  if(iJet->energy()== 0) return false;
347  if(iJet->neutralHadronEnergy()/iJet->energy() > 0.99) return false;
348  if(iJet->neutralEmEnergy()/iJet->energy() > 0.99) return false;
349  if(iJet->nConstituents() < 2) return false;
350  if(iJet->chargedHadronEnergy()/iJet->energy() <= 0 && std::abs(iJet->eta()) < 2.4 ) return false;
351  if(iJet->chargedEmEnergy()/iJet->energy() > 0.99 && std::abs(iJet->eta()) < 2.4 ) return false;
352  if(iJet->chargedMultiplicity() < 1 && std::abs(iJet->eta()) < 2.4 ) return false;
353  return true;
354 }
355 
358  const reco::Vertex* hardScatterVertex,double iDRMax)
359 {
360  double iDR2Max = iDRMax*iDRMax;
361  reco::Candidate::LorentzVector lVis(0,0,0,0);
362  for ( reco::CandidateView::const_iterator pfCandidate = pfCandidates.begin();
363  pfCandidate != pfCandidates.end(); ++pfCandidate ) {
364  if(deltaR2(iCand->p4(),pfCandidate->p4()) > iDR2Max) continue;
365  double dZ = -999.; // PH: If no vertex is reconstructed in the event
366  // or PFCandidate has no track, set dZ to -999
367  if ( hardScatterVertex ) {
368  const reco::PFCandidate* pfc = dynamic_cast<const reco::PFCandidate* >( (&(*pfCandidate)) );
369  if( pfc != nullptr ) { //PF candidate for RECO and PAT levels
370  if ( pfc->trackRef().isNonnull() ) dZ = std::abs(pfc->trackRef()->dz(hardScatterVertex->position()));
371  else if ( pfc->gsfTrackRef().isNonnull() ) dZ = std::abs(pfc->gsfTrackRef()->dz(hardScatterVertex->position()));
372  }
373  else { //if not, then packedCandidate for miniAOD level
374  const pat::PackedCandidate* pfc = dynamic_cast<const pat::PackedCandidate* >( &(*pfCandidate) );
375  dZ = std::abs( pfc->dz( hardScatterVertex->position() ) );
376  }
377  }
378  if(std::abs(dZ) > 0.1) continue;
379  lVis += pfCandidate->p4();
380  }
381  return lVis.pt()/iCand->pt();
382 }
383 
385 
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:39
virtual float dz(size_t ipv=0) const
dz with respect to the PV[ipv]
PFMETAlgorithmMVA mvaMEtAlgo_
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
const double dR2Max
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
void setDZ(float dZ)
Definition: PUSubMETData.h:51
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:253
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:74
double eta() const final
momentum pseudorapidity
PFMETProducerMVA(const edm::ParameterSet &)
void produce(edm::Event &, const edm::EventSetup &) override
float chargedEmEnergy() const
chargedEmEnergy
Definition: PFJet.h:142
std::vector< edm::EDGetTokenT< reco::CandidateView > > srcLeptons_
PFSpecificAlgo pfMEtSpecificAlgo_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual bool isMuon() const =0
std::vector< reco::Vertex::Point > computeVertexInfo(const reco::VertexCollection &)
virtual bool isPhoton() const =0
double px() const final
x coordinate of momentum vector
void setSignificanceMatrix(const reco::METCovMatrix &matrix)
Definition: MET.cc:157
void setP4(const reco::Candidate::LorentzVector p4)
Definition: PUSubMETData.h:50
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:63
void setMvaVal(float mva)
Definition: PUSubMETData.h:60
size_type size() const
double pt() const final
transverse momentum
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
double correction(const LorentzVector &fJet) const
get correction using Jet information only
Definition: JetCorrector.h:49
reco::Candidate::LorentzVector getMEt() const
reco::CandidatePtrVector signalCands() const
int chargedMultiplicity() const
chargedMultiplicity
Definition: PFJet.h:155
const Point & position() const
position
Definition: Vertex.h:109
const std::vector< reco::PFCandidatePtr > & signalPFCands() const
PFCandidates in signal region.
Definition: PFTau.cc:78
Jets made from PFObjects.
Definition: PFJet.h:21
void setInput(const std::vector< reco::PUSubMETCandInfo > &, const std::vector< reco::PUSubMETCandInfo > &, const std::vector< reco::PUSubMETCandInfo > &, const std::vector< reco::Vertex::Point > &)
double chargedFracInCone(const reco::Candidate *iCand, const reco::CandidateView &pfCandidates, const reco::Vertex *hardScatterVertex, double iDRMax=0.2)
float neutralEmEnergy() const
neutralEmEnergy
Definition: PFJet.h:150
const reco::METCovMatrix & getMEtCov() const
edm::EDGetTokenT< edm::ValueMap< float > > srcJetIds_
reco::TrackRef trackRef() const
Definition: PFCandidate.cc:442
bool istau(const reco::Candidate *iCand)
int iEvent
Definition: GenABIO.cc:230
const_iterator begin() const
edm::EDGetTokenT< edm::View< reco::Candidate > > srcPFCandidatesView_
CommonMETData run(const edm::View< reco::Candidate > &candidates, double globalThreshold=0.0)
Definition: METAlgo.cc:16
bool passPFLooseId(const reco::PFJet *iJet)
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
edm::EDGetTokenT< reco::PFJetCollection > srcCorrJets_
std::vector< edm::InputTag > vInputTag
double p4[4]
Definition: TauolaWrapper.h:92
double energy() const final
energy
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< reco::VertexCollection > srcVertices_
const double dPtMatch
std::vector< reco::PUSubMETCandInfo > computeLeptonInfo(const std::vector< edm::EDGetTokenT< reco::CandidateView > > &srcLeptons_, const reco::CandidateView &pfCandidates, const reco::Vertex *hardScatterVertex, int &lId, bool &lHasPhotons, edm::Event &iEvent)
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
virtual bool isElectron() const =0
void initialize(const edm::EventSetup &)
MET made from Particle Flow Candidates.
Analysis-level tau class.
Definition: Tau.h:55
void setHasPhotons(bool hasPhotons)
SpecificPFMETData run(const edm::View< reco::Candidate > &pfCands)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
void setChargedEnFrac(float chEnF)
Definition: PUSubMETData.h:61
edm::EDGetTokenT< reco::JetCorrector > mJetCorrector_
virtual double pt() const =0
transverse momentum
std::vector< reco::PUSubMETCandInfo > computeJetInfo(const reco::PFJetCollection &, const edm::Handle< reco::PFJetCollection > &, const edm::ValueMap< float > &, const reco::VertexCollection &, const reco::Vertex *, const reco::JetCorrector &iCorr, edm::Event &iEvent, const edm::EventSetup &iSetup, std::vector< reco::PUSubMETCandInfo > &iLeptons, std::vector< reco::PUSubMETCandInfo > &iCands)
double py() const final
y coordinate of momentum vector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
edm::EventID id() const
Definition: EventBase.h:60
const double dR2Min
std::vector< PFJet > PFJetCollection
collection of PFJet objects
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
fixed size matrix
virtual int nConstituents() const
of constituents
Definition: Jet.h:65
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
reco::GsfTrackRef gsfTrackRef() const
Definition: PFCandidate.cc:480
float chargedEnFrac() const
Definition: PUSubMETData.h:45
double chargedEnFrac(const reco::Candidate *iCand, const reco::CandidateView &pfCandidates, const reco::Vertex *hardScatterVertex)
float neutralHadronEnergy() const
neutralHadronEnergy
Definition: PFJet.h:102
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
const_iterator end() const
edm::EDGetTokenT< reco::PFJetCollection > srcUncorrJets_
double phi() const final
momentum azimuthal angle
void setP4(const LorentzVector &p4) final
set 4-momentum
float chargedHadronEnergy() const
chargedHadronEnergy
Definition: PFJet.h:98
def move(src, dest)
Definition: eostools.py:510
const reco::Candidate::LorentzVector & p4() const
Definition: PUSubMETData.h:34
std::vector< reco::PUSubMETCandInfo > computePFCandidateInfo(const reco::CandidateView &, const reco::Vertex *)