CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
HLTJetCollectionsForBoostedLeptonPlusJets< jetType > Class Template Reference

#include <HLTJetCollectionsForBoostedLeptonPlusJets.h>

Inheritance diagram for HLTJetCollectionsForBoostedLeptonPlusJets< jetType >:
edm::stream::EDProducer<>

Public Member Functions

 HLTJetCollectionsForBoostedLeptonPlusJets (const edm::ParameterSet &)
 
 ~HLTJetCollectionsForBoostedLeptonPlusJets () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

void produce (edm::Event &, const edm::EventSetup &) override
 

Private Attributes

edm::InputTag hltLeptonTag
 
edm::EDGetTokenT< std::vector< jetType > > m_theJetToken
 
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefsm_theLeptonToken
 
double minDeltaR2_
 
edm::InputTag sourceJetTag
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

template<typename jetType>
class HLTJetCollectionsForBoostedLeptonPlusJets< jetType >

This class is an EDProducer implementing an HLT trigger for lepton and jet objects, cutting on variables relating to the jet 4-momentum representation. The producer checks for overlaps between leptons and jets and if a combination of one lepton + jets cleaned against this leptons satisfy the cuts. These jets are then added to a cleaned jet collection which is put into the event.

Author
Lukasz Kreczko

Definition at line 29 of file HLTJetCollectionsForBoostedLeptonPlusJets.h.

Constructor & Destructor Documentation

◆ HLTJetCollectionsForBoostedLeptonPlusJets()

template<typename jetType >
HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::HLTJetCollectionsForBoostedLeptonPlusJets ( const edm::ParameterSet iConfig)
explicit

Definition at line 23 of file HLTJetCollectionsForBoostedLeptonPlusJets.cc.

References HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::hltLeptonTag, HLT_2024v10_cff::jetType, HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::m_theJetToken, HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::m_theLeptonToken, and HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::sourceJetTag.

25  : hltLeptonTag(iConfig.getParameter<edm::InputTag>("HltLeptonTag")),
26  sourceJetTag(iConfig.getParameter<edm::InputTag>("SourceJetTag")),
27  // minimum delta-R^2 threshold with sign
28  minDeltaR2_(iConfig.getParameter<double>("minDeltaR") * std::abs(iConfig.getParameter<double>("minDeltaR"))) {
29  using namespace edm;
30  using namespace std;
31 
32  typedef vector<RefVector<vector<jetType>, jetType, refhelper::FindUsingAdvance<vector<jetType>, jetType>>>
33  JetCollectionVector;
34  typedef vector<jetType> JetCollection;
35 
36  m_theLeptonToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltLeptonTag);
37  m_theJetToken = consumes<std::vector<jetType>>(sourceJetTag);
38  produces<JetCollectionVector>();
39  produces<JetCollection>();
40 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::vector< Jet > JetCollection
Definition: Jet.h:53
edm::EDGetTokenT< std::vector< jetType > > m_theJetToken
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_theLeptonToken
HLT enums.

◆ ~HLTJetCollectionsForBoostedLeptonPlusJets()

template<typename jetType >
HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::~HLTJetCollectionsForBoostedLeptonPlusJets ( )
override

Definition at line 43 of file HLTJetCollectionsForBoostedLeptonPlusJets.cc.

43  {
44  // do anything here that needs to be done at desctruction time
45  // (e.g. close files, deallocate resources etc.)
46 }

Member Function Documentation

◆ fillDescriptions()

template<typename jetType >
void HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 49 of file HLTJetCollectionsForBoostedLeptonPlusJets.cc.

References edm::ConfigurationDescriptions::add(), defaultModuleLabel(), submitPVResolutionJobs::desc, and ProducerED_cfi::InputTag.

49  {
51  desc.add<edm::InputTag>("HltLeptonTag", edm::InputTag("triggerFilterObjectWithRefs"));
52  //(2)
53  desc.add<edm::InputTag>("SourceJetTag", edm::InputTag("jetCollection"));
54  //(2)
55  desc.add<double>("minDeltaR", 0.5);
57 }
std::string defaultModuleLabel()
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ produce()

template<typename jetType >
void HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 66 of file HLTJetCollectionsForBoostedLeptonPlusJets.cc.

References reco::deltaR2(), MillePedeFileConverter_cfg::e, spr::find(), trigger::TriggerRefsCollections::getObjects(), iEvent, metsig::jet, PDWG_EXODelayedJetMET_cff::jets, HLT_2024v10_cff::jetType, HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::m_theJetToken, HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::m_theLeptonToken, HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::minDeltaR2_, eostools::move(), jetUpdater_cfi::sort, trigger::TriggerCluster, trigger::TriggerElectron, trigger::TriggerMuon, trigger::TriggerPhoton, and trackerHitRTTI::vector.

66  {
67  using namespace edm;
68  using namespace std;
69  //(3)
70  using namespace reco;
71  //(3)
72  typedef vector<RefVector<vector<jetType>, jetType, refhelper::FindUsingAdvance<vector<jetType>, jetType>>>
73  JetCollectionVector;
74  typedef vector<jetType> JetCollection;
77 
79  iEvent.getByToken(m_theLeptonToken, PrevFilterOutput);
80 
81  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
82  // Electrons can be stored as objects of types TriggerCluster, TriggerElectron, or TriggerPhoton
83  vector<reco::RecoChargedCandidateRef> muonCands;
84  PrevFilterOutput->getObjects(trigger::TriggerMuon, muonCands);
85 
86  vector<Ref<reco::ElectronCollection>> eleCands;
87  PrevFilterOutput->getObjects(trigger::TriggerElectron, eleCands);
88 
89  trigger::VRphoton photonCands;
90  PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
91 
92  vector<Ref<reco::RecoEcalCandidateCollection>> clusCands;
93  PrevFilterOutput->getObjects(trigger::TriggerCluster, clusCands);
94 
95  Handle<JetCollection> theJetCollectionHandle;
96  iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
97 
98  typename JetCollection::const_iterator jet;
99 
100  unique_ptr<JetCollection> allSelections(new JetCollection);
101  unique_ptr<JetCollectionVector> product(new JetCollectionVector);
102 
103  std::vector<size_t> usedCands;
104 
105  if (!muonCands.empty()) { // muons
106  for (jet = theJetCollectionHandle->begin(); jet != theJetCollectionHandle->end(); jet++) {
107  //const jetType* referenceJet = &*jet;
108  jetType cleanedJet = *jet; //copy original jet
109  for (size_t candNr = 0; candNr < muonCands.size(); candNr++) {
110  if (std::find(usedCands.begin(), usedCands.end(), candNr) != usedCands.end())
111  continue;
112  if (reco::deltaR2((*muonCands[candNr]), cleanedJet) <= minDeltaR2_) {
113  std::vector<edm::Ptr<reco::PFCandidate>> pfConstituents = cleanedJet.getPFConstituents();
114  for (std::vector<edm::Ptr<reco::PFCandidate>>::const_iterator i_candidate = pfConstituents.begin();
115  i_candidate != pfConstituents.end();
116  ++i_candidate) {
117  if (reco::deltaR2((*muonCands[candNr]), (**i_candidate)) < 1e-6) {
118  cleanedJet.setP4(cleanedJet.p4() - muonCands[candNr]->p4());
119  usedCands.push_back(candNr);
120  break;
121  } //if constituent matched
122  } //for constituents
123  } //if dR<min
124  } //for cands
125  allSelections->push_back(cleanedJet);
126  } //for jets
127  } //if cands
128 
129  if (!eleCands.empty()) { // electrons
130  for (jet = theJetCollectionHandle->begin(); jet != theJetCollectionHandle->end(); jet++) {
131  //const jetType* referenceJet = &*jet;
132  jetType cleanedJet = *jet; //copy original jet
133  for (size_t candNr = 0; candNr < eleCands.size(); candNr++) {
134  if (std::find(usedCands.begin(), usedCands.end(), candNr) != usedCands.end())
135  continue;
136  if (reco::deltaR2((*eleCands[candNr]), cleanedJet) <= minDeltaR2_) {
137  std::vector<edm::Ptr<reco::PFCandidate>> pfConstituents = cleanedJet.getPFConstituents();
138  for (std::vector<edm::Ptr<reco::PFCandidate>>::const_iterator i_candidate = pfConstituents.begin();
139  i_candidate != pfConstituents.end();
140  ++i_candidate) {
141  if (reco::deltaR2((*eleCands[candNr]), (**i_candidate)) < 1e-6) {
142  cleanedJet.setP4(cleanedJet.p4() - eleCands[candNr]->p4());
143  usedCands.push_back(candNr);
144  break;
145  } //if constituent matched
146  } //for constituents
147  } //if dR<min
148  } //for cands
149  allSelections->push_back(cleanedJet);
150  } //for jets
151  } //if cands
152 
153  if (!photonCands.empty()) { // photons
154  for (jet = theJetCollectionHandle->begin(); jet != theJetCollectionHandle->end(); jet++) {
155  //const jetType* referenceJet = &*jet;
156  jetType cleanedJet = *jet; //copy original jet
157  for (size_t candNr = 0; candNr < photonCands.size(); candNr++) {
158  if (std::find(usedCands.begin(), usedCands.end(), candNr) != usedCands.end())
159  continue;
160  if (reco::deltaR2((*photonCands[candNr]), cleanedJet) <= minDeltaR2_) {
161  std::vector<edm::Ptr<reco::PFCandidate>> pfConstituents = cleanedJet.getPFConstituents();
162  for (std::vector<edm::Ptr<reco::PFCandidate>>::const_iterator i_candidate = pfConstituents.begin();
163  i_candidate != pfConstituents.end();
164  ++i_candidate) {
165  if (reco::deltaR2((*photonCands[candNr]), (**i_candidate)) < 1e-6) {
166  cleanedJet.setP4(cleanedJet.p4() - photonCands[candNr]->p4());
167  usedCands.push_back(candNr);
168  break;
169  } //if constituent matched
170  } //for constituents
171  } //if dR<min
172  } //for cands
173  allSelections->push_back(cleanedJet);
174  } //for jets
175  } //if cands
176 
177  if (!clusCands.empty()) { // trigger clusters
178  for (jet = theJetCollectionHandle->begin(); jet != theJetCollectionHandle->end(); jet++) {
179  //const jetType* referenceJet = &*jet;
180  jetType cleanedJet = *jet; //copy original jet
181  for (size_t candNr = 0; candNr < clusCands.size(); candNr++) {
182  if (std::find(usedCands.begin(), usedCands.end(), candNr) != usedCands.end())
183  continue;
184  if (reco::deltaR2((*clusCands[candNr]), cleanedJet) <= minDeltaR2_) {
185  std::vector<edm::Ptr<reco::PFCandidate>> pfConstituents = cleanedJet.getPFConstituents();
186  for (std::vector<edm::Ptr<reco::PFCandidate>>::const_iterator i_candidate = pfConstituents.begin();
187  i_candidate != pfConstituents.end();
188  ++i_candidate) {
189  if (reco::deltaR2((*clusCands[candNr]), (**i_candidate)) < 1e-6) {
190  cleanedJet.setP4(cleanedJet.p4() - clusCands[candNr]->p4());
191  usedCands.push_back(candNr);
192  break;
193  } //if constituent matched
194  } //for constituents
195  } //if dR<min
196  } //for cands
197  allSelections->push_back(cleanedJet);
198  } //for jets
199  } //if cands
200 
202  // reorder cleaned jets
203  std::sort(allSelections->begin(), allSelections->end(), compJets);
204  edm::OrphanHandle<JetCollection> cleanedJetHandle = iEvent.put(std::move(allSelections));
205 
206  JetCollection const& jets = *cleanedJetHandle;
207 
208  JetRefVector cleanedJetRefs;
209  cleanedJetRefs.reserve(jets.size());
210  for (unsigned iJet = 0; iJet < jets.size(); ++iJet) {
211  cleanedJetRefs.push_back(JetRef(cleanedJetHandle, iJet));
212  }
213 
214  product->emplace_back(cleanedJetRefs);
215  iEvent.put(std::move(product));
216 
217  return;
218 }
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
std::vector< Jet > JetCollection
Definition: Jet.h:53
edm::Ref< JetBxCollection > JetRef
Definition: Jet.h:12
edm::EDGetTokenT< std::vector< jetType > > m_theJetToken
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_theLeptonToken
edm::RefVector< JetBxCollection > JetRefVector
Definition: Jet.h:13
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
fixed size matrix
HLT enums.
std::vector< reco::RecoEcalCandidateRef > VRphoton
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ hltLeptonTag

template<typename jetType >
edm::InputTag HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::hltLeptonTag
private

◆ m_theJetToken

template<typename jetType >
edm::EDGetTokenT<std::vector<jetType> > HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::m_theJetToken
private

◆ m_theLeptonToken

template<typename jetType >
edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::m_theLeptonToken
private

◆ minDeltaR2_

template<typename jetType >
double HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::minDeltaR2_
private

◆ sourceJetTag

template<typename jetType >
edm::InputTag HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::sourceJetTag
private