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 
4 #include <vector>
5 
10 
13 
14 
15 
29 template <typename JetType>
31 {
32 public:
33  typedef std::vector<JetType> JetCollection;
36 
37  typedef std::vector<edm::RefVector<JetCollection, JetType,
39  //^ This is the type expected by HLTJetCollectionsFilter
40 
41 private:
48  class EtaPhiE
49  {
50  public:
52  EtaPhiE(double eta, double phi, double e);
53 
54  public:
56  double eta() const;
57 
59  double phi() const;
60 
62  double e() const;
63 
65  bool operator<(EtaPhiE const &rhs) const;
66 
67  private:
69  double etaValue, phiValue;
70 
72  double eValue;
73  };
74 
75 public:
78 
79 public:
81  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
82 
84  virtual void produce(edm::Event &iEvent, edm::EventSetup const &iSetup);
85 
86 private:
89 
92 
94  double minDeltaR2;
95 
102  unsigned numLeptons;
103 };
104 
105 
106 
107 // Implementation
108 
110 
115 
117 
118 
119 template <typename JetType>
121  etaValue(eta), phiValue(phi),
122  eValue(e)
123 {}
124 
125 
126 template <typename JetType>
128 {
129  return etaValue;
130 }
131 
132 
133 template <typename JetType>
135 {
136  return phiValue;
137 }
138 
139 
140 template <typename JetType>
142 {
143  return eValue;
144 }
145 
146 
147 template <typename JetType>
149 {
150  return (eValue < rhs.eValue);
151 }
152 
153 
154 template <typename JetType>
156  edm::ParameterSet const &iConfig):
157  minDeltaR2(std::pow(iConfig.getParameter<double>("minDeltaR"), 2)),
158  numLeptons(iConfig.getParameter<unsigned>("numLeptons"))
159 {
160  leptonToken = consumes<trigger::TriggerFilterObjectWithRefs>(
161  iConfig.getParameter<edm::InputTag>("leptons"));
162  jetToken = consumes<std::vector<JetType>>(iConfig.getParameter<edm::InputTag>("jets"));
163 
164  produces<JetCollectionVector>();
165 }
166 
167 
168 template <typename JetType>
170  edm::ConfigurationDescriptions &descriptions)
171 {
173 
174  desc.add<edm::InputTag>("leptons", edm::InputTag("triggerFilterObjectWithRefs"))->
175  setComment("A collection of leptons that pass an HLT filter");
176  desc.add<edm::InputTag>("jets", edm::InputTag("jetCollection"))->
177  setComment("A collection of jets");
178  desc.add<double>("minDeltaR", 0.3)->
179  setComment("Minimal allowed angular separation between a jet and a lepton");
180  desc.add<unsigned>("numLeptons", 1)->
181  setComment("Number of leading leptons against which the jets are cleaned");
182 
184 }
185 
186 
187 template <typename JetType>
189  edm::EventSetup const &iSetup)
190 {
191  // Read results of the lepton filter
193  iEvent.getByToken(leptonToken, filterOutput);
194 
195  // Momenta of the leptons that passed the filter will be pushed into a vector
196  std::vector<EtaPhiE> leptonMomenta;
197 
198 
199  // First, assume these are muons and try to store their momenta
201  filterOutput->getObjects(trigger::TriggerMuon, muons);
202 
203  for (auto const &muRef: muons) // the collection might be empty
204  leptonMomenta.emplace_back(muRef->eta(), muRef->phi(), muRef->energy());
205 
206 
207  // Then get the momenta as if these are electrons. Electrons are tricky because they can be
208  //stored with three types of trigger objects: TriggerElectron, TriggerPhoton, and
209  //TriggerCluster. Try them one by one.
211  filterOutput->getObjects(trigger::TriggerElectron, electrons);
212 
213  for (auto const &eRef: electrons) // the collection might be empty
214  {
215  auto const &sc = eRef->superCluster();
216  leptonMomenta.emplace_back(sc->eta(), sc->phi(), sc->energy());
217  }
218 
220  filterOutput->getObjects(trigger::TriggerPhoton, photons);
221 
222  for (auto const &eRef: photons) // the collection might be empty
223  {
224  auto const &sc = eRef->superCluster();
225  leptonMomenta.emplace_back(sc->eta(), sc->phi(), sc->energy());
226  }
227 
229  filterOutput->getObjects(trigger::TriggerCluster, clusters);
230 
231  for (auto const &eRef: clusters) // the collection might be empty
232  {
233  auto const &sc = eRef->superCluster();
234  leptonMomenta.emplace_back(sc->eta(), sc->phi(), sc->energy());
235  }
236 
237 
238  // Make sure the momenta are sorted
239  std::sort(leptonMomenta.rbegin(), leptonMomenta.rend());
240 
241 
242  // Read the source collection of jets
243  edm::Handle<JetCollection> jetHandle;
244  iEvent.getByToken(jetToken, jetHandle);
245  JetCollection const &jets = *jetHandle;
246 
247 
248  // Put references to jets that are not matched to leptons into a dedicated collection
249  JetRefVector cleanedJetRefs;
250  unsigned const numLeptonsToLoop = std::min<unsigned>(leptonMomenta.size(), numLeptons);
251 
252  for (unsigned iJet = 0; iJet < jets.size(); ++iJet)
253  {
254  bool overlap = false;
255 
256  for (unsigned iLepton = 0; iLepton < numLeptonsToLoop; ++iLepton)
257  if (reco::deltaR2(leptonMomenta.at(iLepton), jets.at(iJet)) < minDeltaR2)
258  {
259  overlap = true;
260  break;
261  }
262 
263  if (not overlap)
264  cleanedJetRefs.push_back(JetRef(jetHandle, iJet));
265  }
266 
267 
268  // Store the collection in the event
269  std::auto_ptr<JetCollectionVector> product(new JetCollectionVector);
270  //^ Have to use the depricated auto_ptr here because this is what edm::Event::put expects
271  product->emplace_back(cleanedJetRefs);
272  iEvent.put(product);
273 }
274 
275 #endif // HLTJetsCleanedFromLeadingLeptons_h
HLTJetsCleanedFromLeadingLeptons(edm::ParameterSet const &iConfig)
Constructor.
T getParameter(std::string const &) const
std::string defaultModuleLabel()
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:462
double phi() const
Returns azimuthal angle.
edm::RefVector< JetCollection > JetRefVector
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.
edm::Ref< JetBxCollection > JetRef
Definition: Jet.h:12
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:121
vector< PseudoJet > jets
double minDeltaR2
A square of the minimal allowed angular separation between a lepton and a jet.
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)
Geom::Phi< T > phi() const
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
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:62
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.