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  produces<reco::PFJetCollection>();
16 }
18 
20  std::unique_ptr<reco::PFJetCollection> cleanedPFJets(new reco::PFJetCollection);
21 
23  iEvent.getByToken(tauSrc_, tauJets);
24 
26  iEvent.getByToken(pfJetSrc_, PFJets);
27 
30 
31  if (PFJets->size() > 1) {
32  for (unsigned int iJet = 0; iJet < PFJets->size(); iJet++) {
33  bool isMatched = false;
34  const reco::PFJet& myPFJet = (*PFJets)[iJet];
35  for (unsigned int iTau = 0; iTau < taus.size(); iTau++) {
36  if (reco::deltaR2(taus[iTau]->p4(), myPFJet.p4()) < matchingR2_) {
37  isMatched = true;
38  break;
39  }
40  }
41  if (isMatched == false)
42  cleanedPFJets->push_back(myPFJet);
43  }
44  }
45  iEvent.put(std::move(cleanedPFJets));
46 }
47 
50  desc.add<edm::InputTag>("PFJetSrc", edm::InputTag("hltAK4PFJetsCorrected"))->setComment("Input collection of PFJets");
51  desc.add<edm::InputTag>("TauSrc", edm::InputTag("hltSinglePFTau20TrackPt1LooseChargedIsolationReg"))
52  ->setComment("Input collection of PFTaus that have passed ID and isolation requirements");
53  desc.add<double>("Min_dR", 0.5)->setComment("Minimum dR outside of which PFJets will be saved");
54  descriptions.setComment(
55  "This module produces a collection of PFJets that are cross-cleaned with respect to PFTaus passing a HLT "
56  "filter.");
57  descriptions.add("PFJetsTauOverlapRemoval", desc);
58 }
59 
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
PFJetsTauOverlapRemoval(const edm::ParameterSet &)
Jets made from PFObjects.
Definition: PFJet.h:20
const LorentzVector & p4() const final
four-momentum Lorentz vector
const edm::EDGetTokenT< reco::PFJetCollection > pfJetSrc_
int iEvent
Definition: GenABIO.cc:224
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
bool isMatched(TrackingRecHit const &hit)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setComment(std::string const &value)
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
void add(std::string const &label, ParameterSetDescription const &psetDescription)
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:511
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)