CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTJetsCleanedFromLeadingLeptons.h
Go to the documentation of this file.
1 #ifndef HLTJetsCleanedFromLeadingLeptons_h
2 #define HLTJetsCleanedFromLeadingLeptons_h
3 
8 
10 
11 #include <vector>
12 
13 
27 template <typename JetType>
29 {
30 public:
31  typedef std::vector<JetType> JetCollection;
34 
35  typedef std::vector<edm::RefVector<JetCollection, JetType,
37  //^ This is the type expected by HLTJetCollectionsFilter
38 
39 private:
46  class EtaPhiE
47  {
48  public:
50  EtaPhiE(double eta, double phi, double e);
51 
52  public:
54  double eta() const;
55 
57  double phi() const;
58 
60  double e() const;
61 
63  bool operator<(EtaPhiE const &rhs) const;
64 
65  private:
67  double etaValue, phiValue;
68 
70  double eValue;
71  };
72 
73 public:
76 
77 public:
79  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
80 
82  virtual void produce(edm::Event &iEvent, edm::EventSetup const &iSetup);
83 
84 private:
87 
90 
92  double minDeltaR2;
93 
100  unsigned numLeptons;
101 };
102 
103 
104 
105 // Implementation
106 
108 
113 
115 
116 
117 template <typename JetType>
119  etaValue(eta), phiValue(phi),
120  eValue(e)
121 {}
122 
123 
124 template <typename JetType>
126 {
127  return etaValue;
128 }
129 
130 
131 template <typename JetType>
133 {
134  return phiValue;
135 }
136 
137 
138 template <typename JetType>
140 {
141  return eValue;
142 }
143 
144 
145 template <typename JetType>
147 {
148  return (eValue < rhs.eValue);
149 }
150 
151 
152 template <typename JetType>
154  edm::ParameterSet const &iConfig):
155  minDeltaR2(std::pow(iConfig.getParameter<double>("minDeltaR"), 2)),
156  numLeptons(iConfig.getParameter<unsigned>("numLeptons"))
157 {
158  leptonToken = consumes<trigger::TriggerFilterObjectWithRefs>(
159  iConfig.getParameter<edm::InputTag>("leptons"));
160  jetToken = consumes<std::vector<JetType>>(iConfig.getParameter<edm::InputTag>("jets"));
161 
162  produces<JetCollectionVector>();
163 }
164 
165 
166 template <typename JetType>
168  edm::ConfigurationDescriptions &descriptions)
169 {
171 
172  desc.add<edm::InputTag>("leptons", edm::InputTag("triggerFilterObjectWithRefs"))->
173  setComment("A collection of leptons that pass an HLT filter");
174  desc.add<edm::InputTag>("jets", edm::InputTag("jetCollection"))->
175  setComment("A collection of jets");
176  desc.add<double>("minDeltaR", 0.3)->
177  setComment("Minimal allowed angular separation between a jet and a lepton");
178  desc.add<unsigned>("numLeptons", 1)->
179  setComment("Number of leading leptons against which the jets are cleaned");
180 
181  descriptions.add(std::string("hlt") +
183 }
184 
185 
186 template <typename JetType>
188  edm::EventSetup const &iSetup)
189 {
190  // Read results of the lepton filter
192  iEvent.getByToken(leptonToken, filterOutput);
193 
194  // Momenta of the leptons that passed the filter will be pushed into a vector
195  std::vector<EtaPhiE> leptonMomenta;
196 
197 
198  // First, assume these are muons and try to store their momenta
200  filterOutput->getObjects(trigger::TriggerMuon, muons);
201 
202  for (auto const &muRef: muons) // the collection might be empty
203  leptonMomenta.emplace_back(muRef->eta(), muRef->phi(), muRef->energy());
204 
205 
206  // Then get the momenta as if these are electrons. Electrons are tricky because they can be
207  //stored with three types of trigger objects: TriggerElectron, TriggerPhoton, and
208  //TriggerCluster. Try them one by one.
210  filterOutput->getObjects(trigger::TriggerElectron, electrons);
211 
212  for (auto const &eRef: electrons) // the collection might be empty
213  {
214  auto const &sc = eRef->superCluster();
215  leptonMomenta.emplace_back(sc->eta(), sc->phi(), sc->energy());
216  }
217 
219  filterOutput->getObjects(trigger::TriggerPhoton, photons);
220 
221  for (auto const &eRef: photons) // the collection might be empty
222  {
223  auto const &sc = eRef->superCluster();
224  leptonMomenta.emplace_back(sc->eta(), sc->phi(), sc->energy());
225  }
226 
227  trigger::VRphoton clusters;
228  filterOutput->getObjects(trigger::TriggerCluster, clusters);
229 
230  for (auto const &eRef: clusters) // the collection might be empty
231  {
232  auto const &sc = eRef->superCluster();
233  leptonMomenta.emplace_back(sc->eta(), sc->phi(), sc->energy());
234  }
235 
236 
237  // Make sure the momenta are sorted
238  std::sort(leptonMomenta.rbegin(), leptonMomenta.rend());
239 
240 
241  // Read the source collection of jets
242  edm::Handle<JetCollection> jetHandle;
243  iEvent.getByToken(jetToken, jetHandle);
244  JetCollection const &jets = *jetHandle;
245 
246 
247  // Put references to jets that are not matched to leptons into a dedicated collection
248  JetRefVector cleanedJetRefs;
249  unsigned const numLeptonsToLoop = std::min<unsigned>(leptonMomenta.size(), numLeptons);
250 
251  for (unsigned iJet = 0; iJet < jets.size(); ++iJet)
252  {
253  bool overlap = false;
254 
255  for (unsigned iLepton = 0; iLepton < numLeptonsToLoop; ++iLepton)
256  if (reco::deltaR2(leptonMomenta.at(iLepton), jets.at(iJet)) < minDeltaR2)
257  {
258  overlap = true;
259  break;
260  }
261 
262  if (not overlap)
263  cleanedJetRefs.push_back(JetRef(jetHandle, iJet));
264  }
265 
266 
267  // Store the collection in the event
268  std::auto_ptr<JetCollectionVector> product(new JetCollectionVector);
269  //^ Have to use the depricated auto_ptr here because this is what edm::Event::put expects
270  product->emplace_back(cleanedJetRefs);
271  iEvent.put(product);
272 }
273 
274 #endif // HLTJetsCleanedFromLeadingLeptons_h
HLTJetsCleanedFromLeadingLeptons(edm::ParameterSet const &iConfig)
Constructor.
T getParameter(std::string const &) const
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > leptonToken
Token to identify a collection of leptons that pass an HLT filter.
unsigned numLeptons
Number of leading leptons against which the jets are cleaned.
bool operator<(EtaPhiE const &rhs) const
A comparison operator to sort a collection of objects of this type.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Describes configuration of the plugin.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
double phi() const
Returns azimuthal angle.
edm::RefVector< JetCollection > JetRefVector
edm::Ref< JetCollection > JetRef
Definition: Jet.h:51
T eta() const
bool overlap(const reco::Muon &muon1, const reco::Muon &muon2, double pullX=1.0, double pullY=1.0, bool checkAdjacentChambers=false)
An auxiliary class to store momentum parametrised in eta, phi, and energy.
virtual void produce(edm::Event &iEvent, edm::EventSetup const &iSetup)
Produces jets cleaned against leptons.
int iEvent
Definition: GenABIO.cc:230
std::vector< edm::RefVector< JetCollection, JetType, edm::refhelper::FindUsingAdvance< JetCollection, JetType > > > JetCollectionVector
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
vector< PseudoJet > jets
double minDeltaR2
A square of the minimal allowed angular separation between a lepton and a jet.
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
EtaPhiE(double eta, double phi, double e)
Constructor.
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Produces a collection of jets cleaned against leading leptons.
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector< reco::RecoChargedCandidateRef > VRmuon
tuple muons
Definition: patZpeak.py:38
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:64
std::vector< reco::ElectronRef > VRelectron
edm::EDGetTokenT< std::vector< JetType > > jetToken
Token to access a collection of jets.
std::vector< reco::RecoEcalCandidateRef > VRphoton
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double etaValue
Pseudorapidity and azimuthal angle.
Definition: DDAxes.h:10