CMS 3D CMS Logo

HLTJetsCleanedFromLeadingLeptons.h
Go to the documentation of this file.
1 #ifndef HLTJetsCleanedFromLeadingLeptons_h
2 #define HLTJetsCleanedFromLeadingLeptons_h
3 
4 #include <vector>
5 
10 
13 
27 template <typename JetType>
29 public:
30  typedef std::vector<JetType> JetCollection;
33 
34  typedef std::vector<edm::RefVector<JetCollection, JetType, edm::refhelper::FindUsingAdvance<JetCollection, JetType>>>
36  //^ This is the type expected by HLTJetCollectionsFilter
37 
38 private:
45  class EtaPhiE {
46  public:
48  EtaPhiE(double eta, double phi, double e);
49 
50  public:
52  double eta() const;
53 
55  double phi() const;
56 
58  double e() const;
59 
61  bool operator<(EtaPhiE const &rhs) const;
62 
63  private:
65  double etaValue, phiValue;
66 
68  double eValue;
69  };
70 
71 public:
74 
75 public:
77  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
78 
80  void produce(edm::Event &iEvent, edm::EventSetup const &iSetup) override;
81 
82 private:
85 
88 
90  double minDeltaR2;
91 
98  unsigned numLeptons;
99 };
100 
101 // Implementation
102 
104 
109 
111 
112 template <typename JetType>
114  : etaValue(eta), phiValue(phi), eValue(e) {}
115 
116 template <typename JetType>
118  return etaValue;
119 }
120 
121 template <typename JetType>
123  return phiValue;
124 }
125 
126 template <typename JetType>
128  return eValue;
129 }
130 
131 template <typename JetType>
133  return (eValue < rhs.eValue);
134 }
135 
136 template <typename JetType>
138  : minDeltaR2(std::pow(iConfig.getParameter<double>("minDeltaR"), 2)),
139  numLeptons(iConfig.getParameter<unsigned>("numLeptons")) {
140  leptonToken = consumes<trigger::TriggerFilterObjectWithRefs>(iConfig.getParameter<edm::InputTag>("leptons"));
141  jetToken = consumes<std::vector<JetType>>(iConfig.getParameter<edm::InputTag>("jets"));
142 
143  produces<JetCollectionVector>();
144 }
145 
146 template <typename JetType>
149 
150  desc.add<edm::InputTag>("leptons", edm::InputTag("triggerFilterObjectWithRefs"))
151  ->setComment("A collection of leptons that pass an HLT filter");
152  desc.add<edm::InputTag>("jets", edm::InputTag("jetCollection"))->setComment("A collection of jets");
153  desc.add<double>("minDeltaR", 0.3)->setComment("Minimal allowed angular separation between a jet and a lepton");
154  desc.add<unsigned>("numLeptons", 1)->setComment("Number of leading leptons against which the jets are cleaned");
155 
157 }
158 
159 template <typename JetType>
161  // Read results of the lepton filter
163  iEvent.getByToken(leptonToken, filterOutput);
164 
165  // Momenta of the leptons that passed the filter will be pushed into a vector
166  std::vector<EtaPhiE> leptonMomenta;
167 
168  // First, assume these are muons and try to store their momenta
170  filterOutput->getObjects(trigger::TriggerMuon, muons);
171 
172  for (auto const &muRef : muons) // the collection might be empty
173  leptonMomenta.emplace_back(muRef->eta(), muRef->phi(), muRef->energy());
174 
175  // Then get the momenta as if these are electrons. Electrons are tricky because they can be
176  //stored with three types of trigger objects: TriggerElectron, TriggerPhoton, and
177  //TriggerCluster. Try them one by one.
180 
181  for (auto const &eRef : electrons) // the collection might be empty
182  {
183  auto const &sc = eRef->superCluster();
184  leptonMomenta.emplace_back(sc->eta(), sc->phi(), sc->energy());
185  }
186 
188  filterOutput->getObjects(trigger::TriggerPhoton, photons);
189 
190  for (auto const &eRef : photons) // the collection might be empty
191  {
192  auto const &sc = eRef->superCluster();
193  leptonMomenta.emplace_back(sc->eta(), sc->phi(), sc->energy());
194  }
195 
198 
199  for (auto const &eRef : clusters) // the collection might be empty
200  {
201  auto const &sc = eRef->superCluster();
202  leptonMomenta.emplace_back(sc->eta(), sc->phi(), sc->energy());
203  }
204 
205  // Make sure the momenta are sorted
206  std::sort(leptonMomenta.rbegin(), leptonMomenta.rend());
207 
208  // Read the source collection of jets
209  edm::Handle<JetCollection> jetHandle;
210  iEvent.getByToken(jetToken, jetHandle);
211  JetCollection const &jets = *jetHandle;
212 
213  // Put references to jets that are not matched to leptons into a dedicated collection
214  JetRefVector cleanedJetRefs;
215  unsigned const numLeptonsToLoop = std::min<unsigned>(leptonMomenta.size(), numLeptons);
216 
217  for (unsigned iJet = 0; iJet < jets.size(); ++iJet) {
218  bool overlap = false;
219 
220  for (unsigned iLepton = 0; iLepton < numLeptonsToLoop; ++iLepton)
221  if (reco::deltaR2(leptonMomenta.at(iLepton), jets.at(iJet)) < minDeltaR2) {
222  overlap = true;
223  break;
224  }
225 
226  if (not overlap)
227  cleanedJetRefs.push_back(JetRef(jetHandle, iJet));
228  }
229 
230  // Store the collection in the event
231  std::unique_ptr<JetCollectionVector> product(new JetCollectionVector);
232  //^ Have to use the depricated unique_ptr here because this is what edm::Event::put expects
233  product->emplace_back(cleanedJetRefs);
234  iEvent.put(std::move(product));
235 }
236 
237 #endif // HLTJetsCleanedFromLeadingLeptons_h
defaultModuleLabel.h
ConfigurationDescriptions.h
HLTJetsCleanedFromLeadingLeptons::EtaPhiE::EtaPhiE
EtaPhiE(double eta, double phi, double e)
Constructor.
Definition: HLTJetsCleanedFromLeadingLeptons.h:113
PDWG_BPHSkim_cff.muons
muons
Definition: PDWG_BPHSkim_cff.py:47
HLTJetsCleanedFromLeadingLeptons::jetToken
edm::EDGetTokenT< std::vector< JetType > > jetToken
Token to access a collection of jets.
Definition: HLTJetsCleanedFromLeadingLeptons.h:87
HLTJetsCleanedFromLeadingLeptons::EtaPhiE::phi
double phi() const
Returns azimuthal angle.
Definition: HLTJetsCleanedFromLeadingLeptons.h:122
HLTJetsCleanedFromLeadingLeptons::EtaPhiE::eValue
double eValue
Energy.
Definition: HLTJetsCleanedFromLeadingLeptons.h:68
Handle.h
HLTJetsCleanedFromLeadingLeptons::EtaPhiE::phiValue
double phiValue
Definition: HLTJetsCleanedFromLeadingLeptons.h:65
trigger::TriggerElectron
Definition: TriggerTypeDefs.h:78
HLTJetsCleanedFromLeadingLeptons::HLTJetsCleanedFromLeadingLeptons
HLTJetsCleanedFromLeadingLeptons(edm::ParameterSet const &iConfig)
Constructor.
Definition: HLTJetsCleanedFromLeadingLeptons.h:137
HLTJetsCleanedFromLeadingLeptons::JetCollection
std::vector< JetType > JetCollection
Definition: HLTJetsCleanedFromLeadingLeptons.h:30
trigger::VRphoton
std::vector< reco::RecoEcalCandidateRef > VRphoton
Definition: TriggerRefsCollections.h:67
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs >
HLTJetsCleanedFromLeadingLeptons::leptonToken
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > leptonToken
Token to identify a collection of leptons that pass an HLT filter.
Definition: HLTJetsCleanedFromLeadingLeptons.h:84
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89301
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
HLTJetsCleanedFromLeadingLeptons
Produces a collection of jets cleaned against leading leptons.
Definition: HLTJetsCleanedFromLeadingLeptons.h:28
EDProducer.h
TriggerFilterObjectWithRefs.h
singleTopDQM_cfi.jets
jets
Definition: singleTopDQM_cfi.py:42
edm::RefVector
Definition: EDProductfwd.h:27
edm::Handle< trigger::TriggerFilterObjectWithRefs >
HLTJetsCleanedFromLeadingLeptons::JetCollectionVector
std::vector< edm::RefVector< JetCollection, JetType, edm::refhelper::FindUsingAdvance< JetCollection, JetType > > > JetCollectionVector
Definition: HLTJetsCleanedFromLeadingLeptons.h:35
HLTJetsCleanedFromLeadingLeptons::EtaPhiE
An auxiliary class to store momentum parametrised in eta, phi, and energy.
Definition: HLTJetsCleanedFromLeadingLeptons.h:45
trigger::TriggerRefsCollections::getObjects
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
Definition: TriggerRefsCollections.h:590
HLTJetsCleanedFromLeadingLeptons::produce
void produce(edm::Event &iEvent, edm::EventSetup const &iSetup) override
Produces jets cleaned against leptons.
Definition: HLTJetsCleanedFromLeadingLeptons.h:160
HLTJetsCleanedFromLeadingLeptons::EtaPhiE::operator<
bool operator<(EtaPhiE const &rhs) const
A comparison operator to sort a collection of objects of this type.
Definition: HLTJetsCleanedFromLeadingLeptons.h:132
edm::Ref
Definition: AssociativeIterator.h:58
HLTJetsCleanedFromLeadingLeptons::EtaPhiE::etaValue
double etaValue
Pseudorapidity and azimuthal angle.
Definition: HLTJetsCleanedFromLeadingLeptons.h:65
HLTJetsCleanedFromLeadingLeptons::EtaPhiE::e
double e() const
Returns energy.
Definition: HLTJetsCleanedFromLeadingLeptons.h:127
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
HLTJetsCleanedFromLeadingLeptons::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Describes configuration of the plugin.
Definition: HLTJetsCleanedFromLeadingLeptons.h:147
PVValHelper::eta
Definition: PVValidationHelpers.h:70
trigger::TriggerMuon
Definition: TriggerTypeDefs.h:79
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
HLTJetsCleanedFromLeadingLeptons::minDeltaR2
double minDeltaR2
A square of the minimal allowed angular separation between a lepton and a jet.
Definition: HLTJetsCleanedFromLeadingLeptons.h:90
bsc_activity_cfg.clusters
clusters
Definition: bsc_activity_cfg.py:36
edm::ParameterSet
Definition: ParameterSet.h:47
defaultModuleLabel
std::string defaultModuleLabel()
Definition: defaultModuleLabel.h:16
Event.h
PVValHelper::phi
Definition: PVValidationHelpers.h:69
deltaR.h
HLTJetsCleanedFromLeadingLeptons::JetRef
edm::Ref< JetCollection > JetRef
Definition: HLTJetsCleanedFromLeadingLeptons.h:31
jetUpdater_cfi.sort
sort
Definition: jetUpdater_cfi.py:29
reco::deltaR2
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
iEvent
int iEvent
Definition: GenABIO.cc:224
muon::overlap
bool overlap(const reco::Muon &muon1, const reco::Muon &muon2, double pullX=1.0, double pullY=1.0, bool checkAdjacentChambers=false)
Definition: MuonSelectors.cc:791
edm::stream::EDProducer
Definition: EDProducer.h:36
BPHMonitor_cfi.photons
photons
Definition: BPHMonitor_cfi.py:91
edm::EventSetup
Definition: EventSetup.h:58
trigger::TriggerCluster
Definition: TriggerTypeDefs.h:88
edm::RefVector::push_back
void push_back(value_type const &ref)
Add a Ref<C, T> to the RefVector.
Definition: RefVector.h:67
l1t::JetRef
edm::Ref< JetBxCollection > JetRef
Definition: Jet.h:12
trigger::VRelectron
std::vector< reco::ElectronRef > VRelectron
Definition: TriggerRefsCollections.h:68
HLTJetsCleanedFromLeadingLeptons::JetRefVector
edm::RefVector< JetCollection > JetRefVector
Definition: HLTJetsCleanedFromLeadingLeptons.h:32
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
Electron.h
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
RecoEcalCandidate.h
pwdgSkimBPark_cfi.electrons
electrons
Definition: pwdgSkimBPark_cfi.py:6
RecoChargedCandidate.h
trigger::VRmuon
std::vector< reco::RecoChargedCandidateRef > VRmuon
Definition: TriggerRefsCollections.h:69
SuperCluster.h
HLTJetsCleanedFromLeadingLeptons::numLeptons
unsigned numLeptons
Number of leading leptons against which the jets are cleaned.
Definition: HLTJetsCleanedFromLeadingLeptons.h:98
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
HLT_FULL_cff.numLeptons
numLeptons
Definition: HLT_FULL_cff.py:26297
HLTJetsCleanedFromLeadingLeptons::EtaPhiE::eta
double eta() const
Returns pseudorapidity.
Definition: HLTJetsCleanedFromLeadingLeptons.h:117
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
trigger::TriggerPhoton
HLT.
Definition: TriggerTypeDefs.h:77
ParameterSet.h
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
unpackData-CaloStage1.jetToken
jetToken
Definition: unpackData-CaloStage1.py:162
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37