CMS 3D CMS Logo

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

#include <HLTJetCollForElePlusJets.h>

Inheritance diagram for HLTJetCollForElePlusJets< T >:
edm::stream::EDProducer<>

Public Member Functions

 HLTJetCollForElePlusJets (const edm::ParameterSet &)
 
 ~HLTJetCollForElePlusJets () override=default
 
- 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 hltElectronTag
 
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefsm_theElectronToken
 
edm::EDGetTokenT< std::vector< T > > m_theJetToken
 
double maxAbsJetEta_
 
double minDeltaEta_
 
double minDeltaR2_
 
double minJetPt_
 
unsigned int minNJets_
 
double minSoftJetPt_
 
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 T>
class HLTJetCollForElePlusJets< T >

This class is an EDProducer implementing an HLT trigger for electron and jet objects, cutting on variables relating to the jet 4-momentum representation. The producer checks for overlaps between electrons and jets and if a combination of one electron + jets cleaned against this electrons 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 30 of file HLTJetCollForElePlusJets.h.

Constructor & Destructor Documentation

◆ HLTJetCollForElePlusJets()

template<typename T >
HLTJetCollForElePlusJets< T >::HLTJetCollForElePlusJets ( const edm::ParameterSet iConfig)
explicit

Definition at line 18 of file HLTJetCollForElePlusJets.cc.

References HLTJetCollForElePlusJets< T >::hltElectronTag, HLTJetCollForElePlusJets< T >::m_theElectronToken, HLTJetCollForElePlusJets< T >::m_theJetToken, and HLTJetCollForElePlusJets< T >::sourceJetTag.

19  : hltElectronTag(iConfig.getParameter<edm::InputTag>("HltElectronTag")),
20  sourceJetTag(iConfig.getParameter<edm::InputTag>("SourceJetTag")),
21  minJetPt_(iConfig.getParameter<double>("MinJetPt")),
22  maxAbsJetEta_(iConfig.getParameter<double>("MaxAbsJetEta")),
23  minNJets_(iConfig.getParameter<unsigned int>("MinNJets")),
24  // minimum delta-R^2 threshold with sign
25  minDeltaR2_(iConfig.getParameter<double>("minDeltaR") * std::abs(iConfig.getParameter<double>("minDeltaR"))),
26  //Only for VBF
27  minSoftJetPt_(iConfig.getParameter<double>("MinSoftJetPt")),
28  minDeltaEta_(iConfig.getParameter<double>("MinDeltaEta")) {
29  typedef std::vector<T> TCollection;
30  m_theElectronToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltElectronTag);
31  m_theJetToken = consumes<TCollection>(sourceJetTag);
32  produces<TCollection>();
33 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_theElectronToken
edm::EDGetTokenT< std::vector< T > > m_theJetToken
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ ~HLTJetCollForElePlusJets()

template<typename T >
HLTJetCollForElePlusJets< T >::~HLTJetCollForElePlusJets ( )
overridedefault

Member Function Documentation

◆ fillDescriptions()

template<typename T >
void HLTJetCollForElePlusJets< T >::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 36 of file HLTJetCollForElePlusJets.cc.

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

36  {
38  desc.add<edm::InputTag>("HltElectronTag", edm::InputTag("triggerFilterObjectWithRefs"));
39  desc.add<edm::InputTag>("SourceJetTag", edm::InputTag("jetCollection"));
40  desc.add<double>("MinJetPt", 30.);
41  desc.add<double>("MaxAbsJetEta", 2.6);
42  desc.add<unsigned int>("MinNJets", 1);
43  desc.add<double>("minDeltaR", 0.5);
44  //Only for VBF
45  desc.add<double>("MinSoftJetPt", 25.);
46  desc.add<double>("MinDeltaEta", -1.);
48 }
std::string defaultModuleLabel()
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ produce()

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

Definition at line 56 of file HLTJetCollForElePlusJets.cc.

References funct::abs(), trigger::TriggerRefsCollections::getObjects(), iEvent, dqmiolumiharvest::j, dqmdumpme::k, eostools::move(), multPhiCorr_741_25nsDY_cfi::px, multPhiCorr_741_25nsDY_cfi::py, jetUpdater_cfi::sort, trigger::TriggerCluster, trigger::TriggerElectron, trigger::TriggerPhoton, and tier0::unique().

56  {
57  using namespace edm;
58  using namespace std;
59 
60  typedef vector<T> TCollection;
61  typedef Ref<TCollection> TRef;
62  typedef edm::RefVector<TCollection> TRefVector;
63  typedef std::vector<edm::RefVector<std::vector<T>, T, edm::refhelper::FindUsingAdvance<std::vector<T>, T>>>
64  TCollectionVector;
65 
67  iEvent.getByToken(m_theElectronToken, PrevFilterOutput);
68 
69  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
70  // Electrons can be stored as objects of types TriggerCluster, TriggerElectron, or TriggerPhoton
71  std::vector<edm::Ref<reco::RecoEcalCandidateCollection>> clusCands;
72  PrevFilterOutput->getObjects(trigger::TriggerCluster, clusCands);
73 
74  std::vector<edm::Ref<reco::ElectronCollection>> eleCands;
75  PrevFilterOutput->getObjects(trigger::TriggerElectron, eleCands);
76 
77  trigger::VRphoton photonCands;
78  PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
79 
80  //prepare the collection of 3-D vector for electron momenta
81  std::vector<TVector3> ElePs;
82 
83  if (!clusCands.empty()) { //try trigger cluster
84  for (auto& clusCand : clusCands) {
85  TVector3 positionVector(clusCand->superCluster()->position().x(),
86  clusCand->superCluster()->position().y(),
87  clusCand->superCluster()->position().z());
88  ElePs.push_back(positionVector);
89  }
90  } else if (!eleCands.empty()) { // try trigger electrons
91  for (auto& eleCand : eleCands) {
92  TVector3 positionVector(eleCand->superCluster()->position().x(),
93  eleCand->superCluster()->position().y(),
94  eleCand->superCluster()->position().z());
95  ElePs.push_back(positionVector);
96  }
97  } else if (!photonCands.empty()) { // try trigger photons
98  for (auto& photonCand : photonCands) {
99  TVector3 positionVector(photonCand->superCluster()->position().x(),
100  photonCand->superCluster()->position().y(),
101  photonCand->superCluster()->position().z());
102  ElePs.push_back(positionVector);
103  }
104  }
105 
106  edm::Handle<TCollection> theJetCollectionHandle;
107  iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
108 
109  const TCollection& theJetCollection = *theJetCollectionHandle;
110 
111  std::unique_ptr<TCollection> theFilteredJetCollection(new TCollection);
112 
113  std::unique_ptr<TCollectionVector> allSelections(new TCollectionVector());
114 
115  bool foundSolution(false);
116 
117  for (auto const& EleP : ElePs) {
118  bool VBFJetPair = false;
119  std::vector<int> store_jet;
120  TRefVector refVector;
121 
122  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
123  TVector3 JetP(theJetCollection[j].px(), theJetCollection[j].py(), theJetCollection[j].pz());
124 
125  double const Deta = EleP.Eta() - JetP.Eta();
126  double const Dphi = EleP.DeltaPhi(JetP);
127  double const DR2 = Deta * Deta + Dphi * Dphi;
128 
129  if (JetP.Pt() > minJetPt_ && std::abs(JetP.Eta()) < maxAbsJetEta_ && DR2 > minDeltaR2_) {
130  store_jet.push_back(j);
131  // The VBF part of the filter
132  if (minDeltaEta_ > 0) {
133  for (unsigned int k = j + 1; k < theJetCollection.size(); k++) {
134  TVector3 SoftJetP(theJetCollection[k].px(), theJetCollection[k].py(), theJetCollection[k].pz());
135 
136  if (std::abs(SoftJetP.Eta() - JetP.Eta()) <= minDeltaEta_) {
137  continue;
138  }
139 
140  double const softDeta = EleP.Eta() - SoftJetP.Eta();
141  double const softDphi = EleP.DeltaPhi(SoftJetP);
142  double const softDR2 = softDeta * softDeta + softDphi * softDphi;
143 
144  if (SoftJetP.Pt() > minSoftJetPt_ && std::abs(SoftJetP.Eta()) < maxAbsJetEta_ && softDR2 > minDeltaR2_) {
145  store_jet.push_back(k);
146  VBFJetPair = true;
147  }
148  }
149  }
150  }
151  }
152 
153  // Now remove duplicates from the jet collection to store
154  std::sort(store_jet.begin(), store_jet.end());
155  store_jet.erase(unique(store_jet.begin(), store_jet.end()), store_jet.end());
156 
157  // Now save the cleaned jets
158  for (int& ijet : store_jet) {
159  //store all selections
160  refVector.push_back(TRef(theJetCollectionHandle, ijet));
161  //store first selection which matches the criteria
162  if (!foundSolution)
163  theFilteredJetCollection->push_back(theJetCollection[ijet]);
164  }
165  //store all selections
166  allSelections->push_back(refVector);
167 
168  if (theFilteredJetCollection->size() >= minNJets_ && minDeltaEta_ < 0)
169  foundSolution = true;
170  else if (VBFJetPair && minDeltaEta_ > 0)
171  foundSolution = true;
172  else if (!foundSolution)
173  theFilteredJetCollection->clear();
174  }
175 
176  iEvent.put(std::move(theFilteredJetCollection));
177 
178  return;
179 }
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_theElectronToken
edm::EDGetTokenT< std::vector< T > > m_theJetToken
int iEvent
Definition: GenABIO.cc:224
def unique(seq, keepstr=True)
Definition: tier0.py:24
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
HLT enums.
long double T
std::vector< reco::RecoEcalCandidateRef > VRphoton
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ hltElectronTag

template<typename T >
edm::InputTag HLTJetCollForElePlusJets< T >::hltElectronTag
private

◆ m_theElectronToken

template<typename T >
edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> HLTJetCollForElePlusJets< T >::m_theElectronToken
private

◆ m_theJetToken

template<typename T >
edm::EDGetTokenT<std::vector<T> > HLTJetCollForElePlusJets< T >::m_theJetToken
private

◆ maxAbsJetEta_

template<typename T >
double HLTJetCollForElePlusJets< T >::maxAbsJetEta_
private

Definition at line 46 of file HLTJetCollForElePlusJets.h.

◆ minDeltaEta_

template<typename T >
double HLTJetCollForElePlusJets< T >::minDeltaEta_
private

Definition at line 50 of file HLTJetCollForElePlusJets.h.

◆ minDeltaR2_

template<typename T >
double HLTJetCollForElePlusJets< T >::minDeltaR2_
private

Definition at line 48 of file HLTJetCollForElePlusJets.h.

◆ minJetPt_

template<typename T >
double HLTJetCollForElePlusJets< T >::minJetPt_
private

Definition at line 45 of file HLTJetCollForElePlusJets.h.

◆ minNJets_

template<typename T >
unsigned int HLTJetCollForElePlusJets< T >::minNJets_
private

Definition at line 47 of file HLTJetCollForElePlusJets.h.

◆ minSoftJetPt_

template<typename T >
double HLTJetCollForElePlusJets< T >::minSoftJetPt_
private

Definition at line 49 of file HLTJetCollForElePlusJets.h.

◆ sourceJetTag

template<typename T >
edm::InputTag HLTJetCollForElePlusJets< T >::sourceJetTag
private