CMS 3D CMS Logo

PFJetsTauOverlapRemoval.cc
Go to the documentation of this file.
2 #include "Math/GenVector/VectorUtil.h"
7 
8 //
9 // class declaration
10 //
12  tauSrc_ ( consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<edm::InputTag>("TauSrc" ) ) ),
13  pfJetSrc_ ( consumes<reco::PFJetCollection>(iConfig.getParameter<edm::InputTag>("PFJetSrc") ) ),
14  matchingR2_ ( iConfig.getParameter<double>("Min_dR")*iConfig.getParameter<double>("Min_dR") )
15 {
16  produces<reco::PFJetCollection>();
17 }
19 
21 {
22  std::unique_ptr<reco::PFJetCollection> cleanedPFJets(new reco::PFJetCollection);
23 
25  iEvent.getByToken(tauSrc_, tauJets);
26 
28  iEvent.getByToken(pfJetSrc_,PFJets);
29 
31  tauJets->getObjects(trigger::TriggerTau,taus);
32 
33  if(PFJets->size() > 1){
34  for(unsigned int iJet = 0; iJet < 2; iJet++){
35  bool isMatched = false;
36  const reco::PFJet & myPFJet = (*PFJets)[iJet];
37  for(unsigned int iTau = 0; iTau < taus.size(); iTau++){
38  if(reco::deltaR2(taus[iTau]->p4(), myPFJet.p4()) < matchingR2_){
39  isMatched = true;
40  break;
41  }
42  }
43  if(isMatched == false) cleanedPFJets->push_back(myPFJet);
44  }
45  if(PFJets->size() > 2) cleanedPFJets->push_back((*PFJets)[2]);
46  }
47  iEvent.put(std::move(cleanedPFJets));
48 }
49 
51 {
53  desc.add<edm::InputTag>("PFJetSrc", edm::InputTag("hltAK4PFJetsCorrected"))->setComment("Input collection of PFJets" );
54  desc.add<edm::InputTag>("TauSrc", edm::InputTag("hltSinglePFTau20TrackPt1LooseChargedIsolationReg"))->setComment("Input collection of PFTaus that have passed ID and isolation requirements");
55  desc.add<double> ("Min_dR",0.5)->setComment("Minimum dR outside of which PFJets will be saved");
56  descriptions.setComment("This module produces a collection of PFJets that are cross-cleaned with respect to PFTaus passing a HLT filter.");
57  descriptions.add ("PFJetsTauOverlapRemoval",desc);
58 }
void setComment(std::string const &value)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
PFJetsTauOverlapRemoval(const edm::ParameterSet &)
Jets made from PFObjects.
Definition: PFJet.h:21
const edm::EDGetTokenT< reco::PFJetCollection > pfJetSrc_
int iEvent
Definition: GenABIO.cc:230
double p4[4]
Definition: TauolaWrapper.h:92
bool isMatched(TrackingRecHit const &hit)
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void setComment(std::string const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
std::vector< PFJet > PFJetCollection
collection of PFJet objects
fixed size matrix
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > tauSrc_
HLT enums.
std::vector< reco::PFTauRef > VRpftau
def move(src, dest)
Definition: eostools.py:510
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override