CMS 3D CMS Logo

HLTDiPFJetPlusTausCandidatePFJetProducer.cc
Go to the documentation of this file.
1 
9 /* Extended for the case of j1, j2, t1 (no t2, i.e. there are only 3 different objects)
10  */
11 
12 #include <algorithm>
13 #include <memory>
14 #include <set>
15 #include <utility>
16 
30 
32 public:
34  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
35  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
36 
37 private:
40  const double extraTauPtCut_;
41  const double mjjMin_, m2jjMin_;
42  const double dRmin_, dRmin2_;
44 };
45 
47  : tauSrc_(consumes(iConfig.getParameter<edm::InputTag>("tauSrc"))),
48  pfJetSrc_(consumes(iConfig.getParameter<edm::InputTag>("pfJetSrc"))),
49  extraTauPtCut_(iConfig.getParameter<double>("extraTauPtCut")),
50  mjjMin_(iConfig.getParameter<double>("mjjMin")),
51  m2jjMin_(mjjMin_ * mjjMin_),
52  dRmin_(iConfig.getParameter<double>("dRmin")),
53  dRmin2_(dRmin_ * dRmin_) {
54  if (dRmin_ <= 0.) {
55  throw cms::Exception("HLTDiPFJetPlusTausCandidatePFJetProducer")
56  << "invalid value for parameter \"dRmin\" (must be > 0): " << dRmin_;
57  }
58  produces<reco::PFJetCollection>();
59 }
60 
63  const edm::EventSetup&) const {
64  auto const& pfJets = iEvent.get(pfJetSrc_);
65 
66  auto cleanedPFJets = std::make_unique<reco::PFJetCollection>();
67  cleanedPFJets->reserve(pfJets.size());
68 
70  iEvent.get(tauSrc_).getObjects(trigger::TriggerTau, taus);
71 
72  std::set<unsigned int> indices;
73 
74  for (unsigned int iJet1 = 0; iJet1 < pfJets.size(); iJet1++) {
75  for (unsigned int iJet2 = iJet1 + 1; iJet2 < pfJets.size(); iJet2++) {
76  bool correctComb = false;
77  const reco::PFJet& myPFJet1 = pfJets[iJet1];
78  const reco::PFJet& myPFJet2 = pfJets[iJet2];
79 
80  if (mjjMin_ >= 0. && (myPFJet1.p4() + myPFJet2.p4()).M2() < m2jjMin_)
81  continue;
82 
83  for (unsigned int iTau1 = 0; iTau1 < taus.size(); iTau1++) {
84  if (reco::deltaR2(taus[iTau1]->p4(), myPFJet1.p4()) < dRmin2_)
85  continue;
86  if (reco::deltaR2(taus[iTau1]->p4(), myPFJet2.p4()) < dRmin2_)
87  continue;
88 
89  if (taus.size() == 1) {
90  if (taus[iTau1]->pt() < extraTauPtCut_)
91  continue;
92 
93  correctComb = true;
94  } else {
95  for (unsigned int iTau2 = iTau1 + 1; iTau2 < taus.size(); iTau2++) {
96  if (taus[iTau1]->pt() < extraTauPtCut_ && taus[iTau2]->pt() < extraTauPtCut_)
97  continue;
98 
99  if (reco::deltaR2(taus[iTau2]->p4(), myPFJet1.p4()) < dRmin2_)
100  continue;
101  if (reco::deltaR2(taus[iTau2]->p4(), myPFJet2.p4()) < dRmin2_)
102  continue;
103 
104  correctComb = true;
105  break;
106  }
107  }
108  if (correctComb)
109  break;
110  }
111 
112  if (correctComb) {
113  indices.insert(iJet1);
114  indices.insert(iJet2);
115  }
116  }
117 
118  for (const auto& i : indices)
119  cleanedPFJets->emplace_back(pfJets[i]);
120  }
121  // sort jets in pt
122  std::sort(cleanedPFJets->begin(), cleanedPFJets->end(), pTComparator_);
123  iEvent.put(std::move(cleanedPFJets));
124 }
125 
128  desc.add<edm::InputTag>("pfJetSrc", edm::InputTag("hltAK4PFJetsCorrected"))->setComment("Input collection of PFJets");
129  desc.add<edm::InputTag>("tauSrc", edm::InputTag("hltSinglePFTau20TrackPt1LooseChargedIsolationReg"))
130  ->setComment("Input collection of PFTaus that have passed ID and isolation requirements");
131  desc.add<double>("extraTauPtCut", 45)->setComment("In case of asymmetric tau pt cuts");
132  desc.add<double>("mjjMin", 500)->setComment("VBF dijet mass condition");
133  desc.add<double>("dRmin", 0.5)->setComment("Minimum dR between PFJets and filtered PFTaus");
134  descriptions.setComment(
135  "This module produces a collection of PFJets that are cross-cleaned with respect to PFTaus passing a HLT "
136  "filter.");
137  descriptions.addWithDefaultLabel(desc);
138 }
139 
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
const edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > tauSrc_
Jets made from PFObjects.
Definition: PFJet.h:20
const LorentzVector & p4() const final
four-momentum Lorentz vector
int iEvent
Definition: GenABIO.cc:224
#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
const edm::EDGetTokenT< reco::PFJetCollection > pfJetSrc_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
HLT enums.
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
std::vector< reco::PFTauRef > VRpftau
def move(src, dest)
Definition: eostools.py:511