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 minDeltaR_
 
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 40 of file HLTJetCollectionsForBoostedLeptonPlusJets.h.

Constructor & Destructor Documentation

◆ HLTJetCollectionsForBoostedLeptonPlusJets()

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

Definition at line 22 of file HLTJetCollectionsForBoostedLeptonPlusJets.cc.

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

24  : hltLeptonTag(iConfig.getParameter<edm::InputTag>("HltLeptonTag")),
25  sourceJetTag(iConfig.getParameter<edm::InputTag>("SourceJetTag")),
26  minDeltaR_(iConfig.getParameter<double>("minDeltaR")) {
27  using namespace edm;
28  using namespace std;
29 
30  typedef vector<RefVector<vector<jetType>, jetType, refhelper::FindUsingAdvance<vector<jetType>, jetType>>>
31  JetCollectionVector;
32  typedef vector<jetType> JetCollection;
33 
34  m_theLeptonToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltLeptonTag);
35  m_theJetToken = consumes<std::vector<jetType>>(sourceJetTag);
36  produces<JetCollectionVector>();
37  produces<JetCollection>();
38 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< Jet > JetCollection
Definition: Jet.h:53
edm::EDGetTokenT< std::vector< jetType > > m_theJetToken
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_theLeptonToken
HLT enums.

◆ ~HLTJetCollectionsForBoostedLeptonPlusJets()

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

Definition at line 41 of file HLTJetCollectionsForBoostedLeptonPlusJets.cc.

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

Member Function Documentation

◆ fillDescriptions()

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

Definition at line 47 of file HLTJetCollectionsForBoostedLeptonPlusJets.cc.

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

47  {
49  desc.add<edm::InputTag>("HltLeptonTag", edm::InputTag("triggerFilterObjectWithRefs"));
50  //(2)
51  desc.add<edm::InputTag>("SourceJetTag", edm::InputTag("jetCollection"));
52  //(2)
53  desc.add<double>("minDeltaR", 0.5);
55 }
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 64 of file HLTJetCollectionsForBoostedLeptonPlusJets.cc.

References PbPb_ZMuSkimMuonDPG_cff::deltaR, spr::find(), trigger::TriggerRefsCollections::getObjects(), iEvent, metsig::jet, PDWG_EXODelayedJetMET_cff::jets, HLT_2022v15_cff::jetType, HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::m_theJetToken, HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::m_theLeptonToken, HLTJetCollectionsForBoostedLeptonPlusJets< jetType >::minDeltaR_, eostools::move(), jetUpdater_cfi::sort, trigger::TriggerCluster, trigger::TriggerElectron, trigger::TriggerMuon, trigger::TriggerPhoton, and trackerHitRTTI::vector.

64  {
65  using namespace edm;
66  using namespace std;
67  //(3)
68  using namespace reco;
69  //(3)
70  typedef vector<RefVector<vector<jetType>, jetType, refhelper::FindUsingAdvance<vector<jetType>, jetType>>>
71  JetCollectionVector;
72  typedef vector<jetType> JetCollection;
75 
77  iEvent.getByToken(m_theLeptonToken, PrevFilterOutput);
78 
79  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
80  // Electrons can be stored as objects of types TriggerCluster, TriggerElectron, or TriggerPhoton
81  vector<reco::RecoChargedCandidateRef> muonCands;
82  PrevFilterOutput->getObjects(trigger::TriggerMuon, muonCands);
83 
84  vector<Ref<reco::ElectronCollection>> eleCands;
85  PrevFilterOutput->getObjects(trigger::TriggerElectron, eleCands);
86 
87  trigger::VRphoton photonCands;
88  PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
89 
90  vector<Ref<reco::RecoEcalCandidateCollection>> clusCands;
91  PrevFilterOutput->getObjects(trigger::TriggerCluster, clusCands);
92 
93  Handle<JetCollection> theJetCollectionHandle;
94  iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
95 
96  typename JetCollection::const_iterator jet;
97 
98  unique_ptr<JetCollection> allSelections(new JetCollection);
99  unique_ptr<JetCollectionVector> product(new JetCollectionVector);
100 
101  std::vector<size_t> usedCands;
102 
103  if (!muonCands.empty()) { // muons
104  for (jet = theJetCollectionHandle->begin(); jet != theJetCollectionHandle->end(); jet++) {
105  //const jetType* referenceJet = &*jet;
106  jetType cleanedJet = *jet; //copy original jet
107  for (size_t candNr = 0; candNr < muonCands.size(); candNr++) {
108  if (std::find(usedCands.begin(), usedCands.end(), candNr) != usedCands.end())
109  continue;
110  if (deltaR((*muonCands[candNr]), cleanedJet) <= minDeltaR_) {
111  std::vector<edm::Ptr<reco::PFCandidate>> pfConstituents = cleanedJet.getPFConstituents();
112  for (std::vector<edm::Ptr<reco::PFCandidate>>::const_iterator i_candidate = pfConstituents.begin();
113  i_candidate != pfConstituents.end();
114  ++i_candidate) {
115  if (deltaR((*muonCands[candNr]), (**i_candidate)) < 0.001) {
116  cleanedJet.setP4(cleanedJet.p4() - muonCands[candNr]->p4());
117  usedCands.push_back(candNr);
118  break;
119  } //if constituent matched
120  } //for constituents
121  } //if dR<min
122  } //for cands
123  allSelections->push_back(cleanedJet);
124  } //for jets
125  } //if cands
126 
127  if (!eleCands.empty()) { // electrons
128  for (jet = theJetCollectionHandle->begin(); jet != theJetCollectionHandle->end(); jet++) {
129  //const jetType* referenceJet = &*jet;
130  jetType cleanedJet = *jet; //copy original jet
131  for (size_t candNr = 0; candNr < eleCands.size(); candNr++) {
132  if (std::find(usedCands.begin(), usedCands.end(), candNr) != usedCands.end())
133  continue;
134  if (deltaR((*eleCands[candNr]), cleanedJet) <= minDeltaR_) {
135  std::vector<edm::Ptr<reco::PFCandidate>> pfConstituents = cleanedJet.getPFConstituents();
136  for (std::vector<edm::Ptr<reco::PFCandidate>>::const_iterator i_candidate = pfConstituents.begin();
137  i_candidate != pfConstituents.end();
138  ++i_candidate) {
139  if (deltaR((*eleCands[candNr]), (**i_candidate)) < 0.001) {
140  cleanedJet.setP4(cleanedJet.p4() - eleCands[candNr]->p4());
141  usedCands.push_back(candNr);
142  break;
143  } //if constituent matched
144  } //for constituents
145  } //if dR<min
146  } //for cands
147  allSelections->push_back(cleanedJet);
148  } //for jets
149  } //if cands
150 
151  if (!photonCands.empty()) { // photons
152  for (jet = theJetCollectionHandle->begin(); jet != theJetCollectionHandle->end(); jet++) {
153  //const jetType* referenceJet = &*jet;
154  jetType cleanedJet = *jet; //copy original jet
155  for (size_t candNr = 0; candNr < photonCands.size(); candNr++) {
156  if (std::find(usedCands.begin(), usedCands.end(), candNr) != usedCands.end())
157  continue;
158  if (deltaR((*photonCands[candNr]), cleanedJet) <= minDeltaR_) {
159  std::vector<edm::Ptr<reco::PFCandidate>> pfConstituents = cleanedJet.getPFConstituents();
160  for (std::vector<edm::Ptr<reco::PFCandidate>>::const_iterator i_candidate = pfConstituents.begin();
161  i_candidate != pfConstituents.end();
162  ++i_candidate) {
163  if (deltaR((*photonCands[candNr]), (**i_candidate)) < 0.001) {
164  cleanedJet.setP4(cleanedJet.p4() - photonCands[candNr]->p4());
165  usedCands.push_back(candNr);
166  break;
167  } //if constituent matched
168  } //for constituents
169  } //if dR<min
170  } //for cands
171  allSelections->push_back(cleanedJet);
172  } //for jets
173  } //if cands
174 
175  if (!clusCands.empty()) { // trigger clusters
176  for (jet = theJetCollectionHandle->begin(); jet != theJetCollectionHandle->end(); jet++) {
177  //const jetType* referenceJet = &*jet;
178  jetType cleanedJet = *jet; //copy original jet
179  for (size_t candNr = 0; candNr < clusCands.size(); candNr++) {
180  if (std::find(usedCands.begin(), usedCands.end(), candNr) != usedCands.end())
181  continue;
182  if (deltaR((*clusCands[candNr]), cleanedJet) <= minDeltaR_) {
183  std::vector<edm::Ptr<reco::PFCandidate>> pfConstituents = cleanedJet.getPFConstituents();
184  for (std::vector<edm::Ptr<reco::PFCandidate>>::const_iterator i_candidate = pfConstituents.begin();
185  i_candidate != pfConstituents.end();
186  ++i_candidate) {
187  if (deltaR((*clusCands[candNr]), (**i_candidate)) < 0.001) {
188  cleanedJet.setP4(cleanedJet.p4() - clusCands[candNr]->p4());
189  usedCands.push_back(candNr);
190  break;
191  } //if constituent matched
192  } //for constituents
193  } //if dR<min
194  } //for cands
195  allSelections->push_back(cleanedJet);
196  } //for jets
197  } //if cands
198 
200  // reorder cleaned jets
201  std::sort(allSelections->begin(), allSelections->end(), compJets);
202  edm::OrphanHandle<JetCollection> cleanedJetHandle = iEvent.put(std::move(allSelections));
203 
204  JetCollection const& jets = *cleanedJetHandle;
205 
206  JetRefVector cleanedJetRefs;
207 
208  for (unsigned iJet = 0; iJet < jets.size(); ++iJet) {
209  cleanedJetRefs.push_back(JetRef(cleanedJetHandle, iJet));
210  }
211 
212  product->emplace_back(cleanedJetRefs);
213  iEvent.put(std::move(product));
214 
215  return;
216 }
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
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

◆ minDeltaR_

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

◆ sourceJetTag

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