CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
TriggerMatcherToHLTDebug Class Reference

#include <MuonAnalysis/MuonAssociators/plugins/TriggerMatcherToHLTDebug.cc>

Inheritance diagram for TriggerMatcherToHLTDebug:
edm::stream::EDProducer<>

Public Member Functions

void produce (edm::Event &event, const edm::EventSetup &eventSetup) override
 
 TriggerMatcherToHLTDebug (const edm::ParameterSet &pset)
 
 ~TriggerMatcherToHLTDebug () override
 Destructor. More...
 
- 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
 

Private Types

typedef edm::AssociationMap
< edm::OneToMany< std::vector
< L2MuonTrajectorySeed >
, std::vector
< L2MuonTrajectorySeed > > > 
SeedMap
 

Private Member Functions

template<typename T >
void storeValueMap (edm::Event &iEvent, const edm::Handle< edm::View< reco::Muon > > &handle, const std::vector< T > &values, const std::string &label) const
 Store extra information in a ValueMap. More...
 

Private Attributes

edm::EDGetTokenT< reco::BeamSpotbeamspotToken_
 
double deltaR_
 
PropagateToMuon l1matcher_
 
edm::EDGetTokenT
< l1extra::L1MuonParticleCollection
l1Token_
 
double max_Dr_L2
 
double max_Dr_L3
 
double max_Dz_L2
 
double max_Dz_L3
 
double max_Eta_L2
 
double max_Eta_L3
 
std::string metname
 
int min_N_L2
 
int min_N_L3
 
int min_Nhits_L2
 
int min_Nhits_L3
 
double min_Pt_L2
 
double min_Pt_L3
 
int minL1Quality_
 
double nsigma_Pt_L2
 
double nsigma_Pt_L3
 
edm::EDGetTokenT< SeedMapseedMapToken_
 
edm::EDGetTokenT< edm::View
< reco::Muon > > 
tagToken_
 
edm::EDGetTokenT
< reco::RecoChargedCandidateCollection
theL2MuonsToken_
 
edm::EDGetTokenT
< L2MuonTrajectorySeedCollection
theL2SeedsToken_
 
edm::EDGetTokenT
< reco::RecoChargedCandidateCollection
theL3MuonsToken_
 
edm::EDGetTokenT
< L3MuonTrajectorySeedCollection
theL3SeedsToken_
 
edm::EDGetTokenT
< reco::TrackCollection
theL3TkTracksToken_
 

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

Description: Matches RECO muons to Trigger ones using HLTDEBUG information. Muon is first matched to L1 using the PropagateToMuon tool from this same package, then all compatible L1s are examined and the corresponding L2 and L3 objects are searched using the references inside those objects.

Definition at line 64 of file TriggerMatcherToHLTDebug.cc.

Member Typedef Documentation

Definition at line 77 of file TriggerMatcherToHLTDebug.cc.

Constructor & Destructor Documentation

TriggerMatcherToHLTDebug::TriggerMatcherToHLTDebug ( const edm::ParameterSet pset)
explicit

Definition at line 131 of file TriggerMatcherToHLTDebug.cc.

References edm::ParameterSet::getParameter(), metname, theL2MuonsToken_, theL2SeedsToken_, theL3MuonsToken_, theL3SeedsToken_, and theL3TkTracksToken_.

132  : tagToken_(consumes<View<reco::Muon> >(pset.getParameter<edm::InputTag>("tags"))),
133  l1Token_(consumes<L1MuonParticleCollection>(pset.getParameter<edm::InputTag>("l1s"))),
134  l1matcher_(pset.getParameter<edm::ParameterSet>("l1matcherConfig"), consumesCollector()),
135  deltaR_(pset.getParameter<double>("deltaR")),
136  minL1Quality_(pset.getParameter<int32_t>("MinL1Quality")),
137  beamspotToken_(consumes<BeamSpot>(pset.getParameter<edm::InputTag>("BeamSpotTag"))),
138  min_N_L2(pset.getParameter<int>("MinN_L2")),
139  max_Eta_L2(pset.getParameter<double>("MaxEta_L2")),
140  min_Nhits_L2(pset.getParameter<int>("MinNhits_L2")),
141  max_Dr_L2(pset.getParameter<double>("MaxDr_L2")),
142  max_Dz_L2(pset.getParameter<double>("MaxDz_L2")),
143  min_Pt_L2(pset.getParameter<double>("MinPt_L2")),
144  nsigma_Pt_L2(pset.getParameter<double>("NSigmaPt_L2")),
145  min_N_L3(pset.getParameter<int>("MinN_L3")),
146  max_Eta_L3(pset.getParameter<double>("MaxEta_L3")),
147  min_Nhits_L3(pset.getParameter<int>("MinNhits_L3")),
148  max_Dr_L3(pset.getParameter<double>("MaxDr_L3")),
149  max_Dz_L3(pset.getParameter<double>("MaxDz_L3")),
150  min_Pt_L3(pset.getParameter<double>("MinPt_L3")),
151  nsigma_Pt_L3(pset.getParameter<double>("NSigmaPt_L3")),
152  seedMapToken_(consumes<SeedMap>(pset.getParameter<edm::InputTag>("SeedMapTag"))) {
153  theL2SeedsToken_ = consumes<L2MuonTrajectorySeedCollection>(pset.getParameter<InputTag>("L2Seeds_Collection"));
154  theL2MuonsToken_ = consumes<RecoChargedCandidateCollection>(pset.getParameter<InputTag>("L2Muons_Collection"));
155  theL3SeedsToken_ = consumes<L3MuonTrajectorySeedCollection>(pset.getParameter<InputTag>("L3Seeds_Collection"));
156  theL3TkTracksToken_ = consumes<TrackCollection>(pset.getParameter<InputTag>("L3TkTracks_Collection"));
157  theL3MuonsToken_ = consumes<RecoChargedCandidateCollection>(pset.getParameter<InputTag>("L3Muons_Collection"));
158 
159  metname = "TriggerMatcherToHLTDebug";
160 
161  produces<edm::ValueMap<int> >("propagatesToM2");
162  produces<edm::ValueMap<int> >("hasL1Particle");
163  produces<edm::ValueMap<int> >("hasL1Filtered");
164  produces<edm::ValueMap<int> >("hasL2Seed");
165  produces<edm::ValueMap<int> >("hasL2Muon");
166  produces<edm::ValueMap<int> >("hasL2MuonFiltered");
167  produces<edm::ValueMap<int> >("hasL3Seed");
168  produces<edm::ValueMap<int> >("hasL3Track");
169  produces<edm::ValueMap<int> >("hasL3Muon");
170  produces<edm::ValueMap<int> >("hasL3MuonFiltered");
171 
172  produces<edm::ValueMap<reco::CandidatePtr> >("l1Candidate");
173  produces<edm::ValueMap<reco::CandidatePtr> >("l2Candidate");
174  produces<edm::ValueMap<reco::CandidatePtr> >("l3Candidate");
175 }
edm::EDGetTokenT< l1extra::L1MuonParticleCollection > l1Token_
edm::EDGetTokenT< L3MuonTrajectorySeedCollection > theL3SeedsToken_
edm::EDGetTokenT< L2MuonTrajectorySeedCollection > theL2SeedsToken_
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > theL2MuonsToken_
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > theL3MuonsToken_
edm::EDGetTokenT< SeedMap > seedMapToken_
edm::EDGetTokenT< edm::View< reco::Muon > > tagToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< reco::BeamSpot > beamspotToken_
edm::EDGetTokenT< reco::TrackCollection > theL3TkTracksToken_
TriggerMatcherToHLTDebug::~TriggerMatcherToHLTDebug ( )
override

Destructor.

Definition at line 178 of file TriggerMatcherToHLTDebug.cc.

178 {}

Member Function Documentation

void TriggerMatcherToHLTDebug::produce ( edm::Event event,
const edm::EventSetup eventSetup 
)
override

Definition at line 181 of file TriggerMatcherToHLTDebug.cc.

References beam_dqm_sourceclient-live_cfg::beamSpot, beamspotToken_, HLT_FULL_cff::deltaR, deltaR_, reco::TrackBase::dxy(), reco::TrackBase::dz(), reco::TrackBase::error(), PV3DBase< T, PVType, FrameType >::eta(), reco::TrackBase::eta(), edmPickEvents::event, PropagateToMuon::extrapolate(), TrajectoryStateOnSurface::globalPosition(), mps_fire::i, edm::HandleBase::id(), edm::Ref< C, T, F >::id(), PropagateToMuon::init(), TrajectoryStateOnSurface::isValid(), edm::HandleBase::isValid(), edm::Ptr< T >::key(), edm::Ref< C, T, F >::key(), l1matcher_, l1Token_, L2Muons_cfi::L2Muons, L3Muons_cfi::L3Muons, LogTrace, max_Dr_L2, max_Dr_L3, max_Dz_L2, max_Dz_L3, max_Eta_L2, max_Eta_L3, metname, min_Nhits_L2, min_Nhits_L3, min_Pt_L2, min_Pt_L3, RPCpg::mu, patZpeak::muons, nsigma_Pt_L2, nsigma_Pt_L3, reco::TrackBase::numberOfValidHits(), reco::TrackBase::parameter(), PV3DBase< T, PVType, FrameType >::phi(), reco::BeamSpot::position(), reco::TrackBase::pt(), L1MuGMTCand::quality(), quality, seedMapToken_, DetachedQuadStep_cff::seeds, edm::RefVector< C, T, F >::size(), tagToken_, theL2MuonsToken_, theL2SeedsToken_, theL3MuonsToken_, theL3SeedsToken_, and theL3TkTracksToken_.

181  {
182  l1matcher_.init(eventSetup);
183 
185  event.getByToken(tagToken_, muons);
186 
188  event.getByToken(l1Token_, L1Muons);
189 
191  event.getByToken(theL2SeedsToken_, L2Seeds);
192 
194  event.getByToken(theL2MuonsToken_, L2Muons);
195 
197  event.getByToken(theL3SeedsToken_, L3Seeds);
198 
200  event.getByToken(theL3TkTracksToken_, L3TkTracks);
201 
203  event.getByToken(theL3MuonsToken_, L3Muons);
204 
205  //beam spot
207  Handle<BeamSpot> recoBeamSpotHandle;
208  event.getByToken(beamspotToken_, recoBeamSpotHandle);
209  beamSpot = *recoBeamSpotHandle;
210 
211  //new for the MAP!!!!
212  edm::Handle<SeedMap> seedMapHandle;
213  event.getByToken(seedMapToken_, seedMapHandle);
214 
215  size_t nmu = muons->size();
216  std::vector<int> propagatesToM2(nmu), hasL1Particle(nmu), hasL1Filtered(nmu);
217  std::vector<int> hasL2Seed(nmu), hasL2Muon(nmu), hasL2MuonFiltered(nmu);
218  std::vector<int> hasL3Seed(nmu), hasL3Track(nmu), hasL3TrackFiltered(nmu), hasL3Muon(nmu), hasL3MuonFiltered(nmu);
219  std::vector<reco::CandidatePtr> l1ptr(nmu), l2ptr(nmu), l3ptr(nmu);
220 
221  for (size_t i = 0; i < nmu; ++i) {
222  const reco::Muon &mu = (*muons)[i];
223 
224  // Propagate to muon station (using the L1 tool)
226  if (!stateAtMB2.isValid())
227  continue;
228  propagatesToM2[i] = 1;
229 
230  double etaTk = stateAtMB2.globalPosition().eta();
231  double phiTk = stateAtMB2.globalPosition().phi();
232  l1extra::L1MuonParticleCollection::const_iterator it;
233  L2MuonTrajectorySeedCollection::const_iterator iSeed;
234  L3MuonTrajectorySeedCollection::const_iterator iSeedL3;
235  RecoChargedCandidateCollection::const_iterator iL2Muon;
236  reco::TrackCollection::const_iterator tktrackL3;
237  RecoChargedCandidateCollection::const_iterator iL3Muon;
238 
239  reco::CandidatePtr thisL1, thisL2, thisL3;
240  for (it = L1Muons->begin(); it != L1Muons->end(); ++it) {
241  const L1MuGMTExtendedCand muonCand = (*it).gmtMuonCand();
242  unsigned int quality = muonCand.quality();
243 
244  double L1phi = (*it).phi();
245  double L1eta = (*it).eta();
246  double L1pt = (*it).pt();
247  double dR = deltaR(etaTk, phiTk, L1eta, L1phi);
248 
249  //CONDIZIONE-> CE NE E' UNA ASSOCIATA?
250  if (dR >= deltaR_)
251  continue;
252  thisL1 = reco::CandidatePtr(L1Muons, it - L1Muons->begin());
253  if (!hasL1Particle[i])
254  l1ptr[i] = thisL1; // if nobody reached L1 before, then we're the best L1 found up to now.
255  hasL1Particle[i]++;
256 
257  if ((quality <= 3) || (L1pt < 7))
258  continue;
259  if (!hasL1Filtered[i])
260  l1ptr[i] = thisL1; // if nobody reached L1 before, then we're the best L1 found up to now.
261  hasL1Filtered[i]++;
262 
263  if (!L2Seeds.isValid())
264  continue;
265  //LOOP SULLA COLLEZIONE DEI SEED
266  for (iSeed = L2Seeds->begin(); iSeed != L2Seeds->end(); ++iSeed) {
267  l1extra::L1MuonParticleRef l1FromSeed = iSeed->l1Particle();
268  if (l1FromSeed.id() != L1Muons.id())
269  throw cms::Exception("CorruptData")
270  << "You're using a different L1 collection than the one used by L2 seeds.\n";
271  if (l1FromSeed.key() != thisL1.key())
272  continue;
273  if (!hasL2Seed[i])
274  l1ptr[i] = thisL1; // if nobody reached here before, we're the best L1
275  hasL2Seed[i]++;
276 
277  if (!L2Muons.isValid())
278  continue;
279  //LOOP SULLA COLLEZIONE L2MUON
280  for (iL2Muon = L2Muons->begin(); iL2Muon != L2Muons->end(); ++iL2Muon) {
281  //MI FACCIO DARE REF E GUARDO SE E' UGUALE AL L2SEED ASSOCIATO
282  //BEFORE THE ASSOCIATION MAP!!!!!
283  //edm::Ref<L2MuonTrajectorySeedCollection> l2seedRef = iL2Muon->track()->seedRef().castTo<edm::Ref<L2MuonTrajectorySeedCollection> >();
284  //l1extra::L1MuonParticleRef l1FromL2 = l2seedRef->l1Particle();
285 
286  //if (l1FromL2.id() != l1FromSeed.id()) throw cms::Exception("CorruptData") << "You're using L2s with a different L1 collection than the one used by L2 seeds.\n";
287  //if (l1FromL2 != l1FromSeed) continue;
288 
289  //AFTER THE ASSOCIATION MAP
291  (*seedMapHandle)[iL2Muon->track()->seedRef().castTo<edm::Ref<L2MuonTrajectorySeedCollection> >()];
292  // bool isTriggered = false;
293  for (size_t jjj = 0; jjj < seeds.size(); jjj++) {
294  if (seeds[jjj]->l1Particle() != l1FromSeed)
295  continue;
296  }
297 
298  thisL2 = reco::CandidatePtr(L2Muons, iL2Muon - L2Muons->begin());
299  if (!hasL2Muon[i]) {
300  l1ptr[i] = thisL1;
301  l2ptr[i] = thisL2;
302  } // if nobody reached here before, we're the best L1 and L2)
303  hasL2Muon[i]++;
304 
305  LogTrace(metname) << "L2MUON TROVATO!" << endl;
306  const reco::Track &L2Track = *iL2Muon->track();
307  double Eta_L2 = L2Track.eta();
308  double Pt_L2 = L2Track.pt();
309  int nValidHits_L2 = L2Track.numberOfValidHits();
310  double BSPos_L2 = L2Track.dxy(beamSpot.position());
311  double dz_L2 = L2Track.dz();
312  double err0_L2 = L2Track.error(0);
313  double abspar0_L2 = fabs(L2Track.parameter(0));
314  double ptLx_L2 = Pt_L2;
315  if (abspar0_L2 > 0)
316  ptLx_L2 += nsigma_Pt_L2 * err0_L2 / abspar0_L2 * Pt_L2;
317 
318  //GUARDO SE L2MUON ASSOCIATO AVREBBE PASSATO IL FILTRO
319  bool passFilter = (((fabs(Eta_L2)) <= max_Eta_L2) && (nValidHits_L2 >= min_Nhits_L2) &&
320  ((fabs(BSPos_L2)) <= max_Dr_L2) && ((fabs(dz_L2)) <= max_Dz_L2) && (ptLx_L2 >= min_Pt_L2));
321  if (!passFilter)
322  continue;
323  if (!hasL2MuonFiltered[i]) {
324  l1ptr[i] = thisL1;
325  l2ptr[i] = thisL2;
326  } // if nobody reached here before, we're the best L1 and L2)
327  hasL2MuonFiltered[i]++;
328 
329  const reco::TrackRef L2FilteredRef = iL2Muon->track();
330 
331  //########L3 PART##############
332  if (!L3Seeds.isValid())
333  continue;
334  for (iSeedL3 = L3Seeds->begin(); iSeedL3 != L3Seeds->end(); ++iSeedL3) {
335  TrackRef staTrack = iSeedL3->l2Track();
336  if (staTrack != L2FilteredRef)
337  continue;
338  if (!hasL3Seed[i]) {
339  l1ptr[i] = thisL1;
340  l2ptr[i] = thisL2;
341  } // if nobody reached here before, we're the best L1 and L2)
342  hasL3Seed[i]++;
343 
344  if (!L3TkTracks.isValid())
345  continue;
346  for (tktrackL3 = L3TkTracks->begin(); tktrackL3 != L3TkTracks->end(); ++tktrackL3) {
348  tktrackL3->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
349  TrackRef staTrack2 = l3seedRef->l2Track();
350 
351  if (staTrack2 != L2FilteredRef)
352  continue;
353  if (!hasL3Track[i]) {
354  l1ptr[i] = thisL1;
355  l2ptr[i] = thisL2;
356  } // if nobody reached here before, we're the best L1 and L2)
357  hasL3Track[i]++;
358 
359  if (!L3Muons.isValid())
360  continue;
361  for (iL3Muon = L3Muons->begin(); iL3Muon != L3Muons->end(); ++iL3Muon) {
363  iL3Muon->track()->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
364  TrackRef staTrack3 = l3seedRef2->l2Track();
365 
366  if (staTrack3 != L2FilteredRef)
367  continue;
368  thisL3 = reco::CandidatePtr(L3Muons, iL3Muon - L3Muons->begin());
369 
370  if (!hasL3Muon[i]) {
371  l1ptr[i] = thisL1;
372  l2ptr[i] = thisL2;
373  l3ptr[i] = thisL3;
374  } // if nobody reached here before, we're the best L1, L2, L3
375  hasL3Muon[i]++;
376 
377  const reco::Track &L3Track = *iL3Muon->track();
378  double Eta_L3 = L3Track.eta();
379  double Pt_L3 = L3Track.pt();
380  int nValidHits_L3 = L3Track.numberOfValidHits();
381  double BSPos_L3 = L3Track.dxy(beamSpot.position());
382  double dz_L3 = L3Track.dz();
383  double err0_L3 = L3Track.error(0);
384  double abspar0_L3 = fabs(L3Track.parameter(0));
385  double ptLx_L3 = Pt_L3;
386 
387  if (abspar0_L3 > 0)
388  ptLx_L3 += nsigma_Pt_L3 * err0_L3 / abspar0_L3 * Pt_L3;
389 
390  if (((fabs(Eta_L3)) <= max_Eta_L3) && (nValidHits_L3 >= min_Nhits_L3) &&
391  ((fabs(BSPos_L3)) <= max_Dr_L3) && ((fabs(dz_L3)) <= max_Dz_L3) && (ptLx_L3 >= min_Pt_L3)) {
392  if (!hasL3MuonFiltered[i]) {
393  l1ptr[i] = thisL1;
394  l2ptr[i] = thisL2;
395  l3ptr[i] = thisL3;
396  } // if nobody reached here before, we're the best L1, L2, L3
397  hasL3MuonFiltered[i]++;
398 
399  } //L3MUON FILTERED ASSOCIATO TROVATO
400  } //L3MUON LOOP
401  } // L3 TRACKS
402  } // L3 SEEDS
403  } //T L2 MUONS
404  } // L2 SEEDS
405  } //L1 MUONS
406  } // RECO MUONS
407  storeValueMap<int>(event, muons, propagatesToM2, "propagatesToM2");
408  storeValueMap<int>(event, muons, hasL1Particle, "hasL1Particle");
409  storeValueMap<int>(event, muons, hasL1Filtered, "hasL1Filtered");
410  storeValueMap<int>(event, muons, hasL2Seed, "hasL2Seed");
411  storeValueMap<int>(event, muons, hasL2Muon, "hasL2Muon");
412  storeValueMap<int>(event, muons, hasL2MuonFiltered, "hasL2MuonFiltered");
413  storeValueMap<int>(event, muons, hasL3Seed, "hasL3Seed");
414  storeValueMap<int>(event, muons, hasL3Track, "hasL3Track");
415  storeValueMap<int>(event, muons, hasL3Muon, "hasL3Muon");
416  storeValueMap<int>(event, muons, hasL3MuonFiltered, "hasL3MuonFiltered");
417  storeValueMap<reco::CandidatePtr>(event, muons, l1ptr, "l1Candidate");
418  storeValueMap<reco::CandidatePtr>(event, muons, l2ptr, "l2Candidate");
419  storeValueMap<reco::CandidatePtr>(event, muons, l3ptr, "l3Candidate");
420 } // METHOD
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
key_type key() const
Definition: Ptr.h:163
ProductID id() const
Definition: HandleBase.cc:29
edm::EDGetTokenT< l1extra::L1MuonParticleCollection > l1Token_
uint32_t const *__restrict__ Quality * quality
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
edm::EDGetTokenT< L3MuonTrajectorySeedCollection > theL3SeedsToken_
edm::EDGetTokenT< L2MuonTrajectorySeedCollection > theL2SeedsToken_
tuple L3Muons
Definition: L3Muons_cfi.py:8
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > theL2MuonsToken_
GlobalPoint globalPosition() const
key_type key() const
Accessor for product key.
Definition: Ref.h:250
tuple L2Muons
Definition: L2Muons_cfi.py:7
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > theL3MuonsToken_
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
#define LogTrace(id)
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
double pt() const
track transverse momentum
Definition: TrackBase.h:637
double error(int i) const
error on specified element
Definition: TrackBase.h:729
edm::EDGetTokenT< SeedMap > seedMapToken_
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:798
const int mu
Definition: Constants.h:22
bool isValid() const
Definition: HandleBase.h:70
double parameter(int i) const
i-th parameter ( i = 0, ... 4 )
Definition: TrackBase.h:723
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:622
edm::EDGetTokenT< edm::View< reco::Muon > > tagToken_
unsigned int quality() const
get quality
Definition: L1MuGMTCand.h:90
edm::Ptr< Candidate > CandidatePtr
persistent reference to an object in a collection of Candidate objects
Definition: CandidateFwd.h:25
void init(const edm::EventSetup &iSetup)
Call this method at the beginning of each event, to initialize geometry, magnetic field and propagato...
T eta() const
Definition: PV3DBase.h:73
tuple muons
Definition: patZpeak.py:39
edm::EDGetTokenT< reco::BeamSpot > beamspotToken_
edm::EDGetTokenT< reco::TrackCollection > theL3TkTracksToken_
size_type size() const
Size of the RefVector.
Definition: RefVector.h:102
const Point & position() const
position
Definition: BeamSpot.h:59
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:608
template<typename T >
void TriggerMatcherToHLTDebug::storeValueMap ( edm::Event iEvent,
const edm::Handle< edm::View< reco::Muon > > &  handle,
const std::vector< T > &  values,
const std::string &  label 
) const
private

Store extra information in a ValueMap.

Definition at line 423 of file TriggerMatcherToHLTDebug.cc.

References edm::helper::Filler< Map >::fill(), patZpeak::handle, edm::helper::Filler< Map >::insert(), eostools::move(), and edm::Event::put().

426  {
427  using namespace edm;
428  using namespace std;
429  unique_ptr<ValueMap<T> > valMap(new ValueMap<T>());
430  typename edm::ValueMap<T>::Filler filler(*valMap);
431  filler.insert(handle, values.begin(), values.end());
432  filler.fill();
433  iEvent.put(std::move(valMap), label);
434 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
char const * label
def move
Definition: eostools.py:511

Member Data Documentation

edm::EDGetTokenT<reco::BeamSpot> TriggerMatcherToHLTDebug::beamspotToken_
private

Definition at line 92 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::deltaR_
private

Definition at line 86 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

PropagateToMuon TriggerMatcherToHLTDebug::l1matcher_
private

Definition at line 81 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

edm::EDGetTokenT<l1extra::L1MuonParticleCollection> TriggerMatcherToHLTDebug::l1Token_
private

Definition at line 80 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Dr_L2
private

Definition at line 96 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Dr_L3
private

Definition at line 105 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Dz_L2
private

Definition at line 97 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Dz_L3
private

Definition at line 106 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Eta_L2
private

Definition at line 94 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Eta_L3
private

Definition at line 103 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

std::string TriggerMatcherToHLTDebug::metname
private

Definition at line 83 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

int TriggerMatcherToHLTDebug::min_N_L2
private

Definition at line 93 of file TriggerMatcherToHLTDebug.cc.

int TriggerMatcherToHLTDebug::min_N_L3
private

Definition at line 102 of file TriggerMatcherToHLTDebug.cc.

int TriggerMatcherToHLTDebug::min_Nhits_L2
private

Definition at line 95 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

int TriggerMatcherToHLTDebug::min_Nhits_L3
private

Definition at line 104 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::min_Pt_L2
private

Definition at line 98 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::min_Pt_L3
private

Definition at line 107 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

int TriggerMatcherToHLTDebug::minL1Quality_
private

Definition at line 89 of file TriggerMatcherToHLTDebug.cc.

double TriggerMatcherToHLTDebug::nsigma_Pt_L2
private

Definition at line 99 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::nsigma_Pt_L3
private

Definition at line 108 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

edm::EDGetTokenT<SeedMap> TriggerMatcherToHLTDebug::seedMapToken_
private

Definition at line 115 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

edm::EDGetTokenT<edm::View<reco::Muon> > TriggerMatcherToHLTDebug::tagToken_
private

Definition at line 79 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

edm::EDGetTokenT<reco::RecoChargedCandidateCollection> TriggerMatcherToHLTDebug::theL2MuonsToken_
private

Definition at line 111 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

edm::EDGetTokenT<L2MuonTrajectorySeedCollection> TriggerMatcherToHLTDebug::theL2SeedsToken_
private

Definition at line 110 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

edm::EDGetTokenT<reco::RecoChargedCandidateCollection> TriggerMatcherToHLTDebug::theL3MuonsToken_
private

Definition at line 114 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

edm::EDGetTokenT<L3MuonTrajectorySeedCollection> TriggerMatcherToHLTDebug::theL3SeedsToken_
private

Definition at line 112 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

edm::EDGetTokenT<reco::TrackCollection> TriggerMatcherToHLTDebug::theL3TkTracksToken_
private

Definition at line 113 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().