CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

void beginRun (const edm::Run &run, const edm::EventSetup &eventSetup) override
 
void produce (edm::Event &event, const edm::EventSetup &eventSetup) override
 
 TriggerMatcherToHLTDebug (const edm::ParameterSet &pset)
 
virtual ~TriggerMatcherToHLTDebug ()
 Destructor. More...
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

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::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

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 72 of file TriggerMatcherToHLTDebug.cc.

Member Typedef Documentation

Definition at line 86 of file TriggerMatcherToHLTDebug.cc.

Constructor & Destructor Documentation

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

Definition at line 141 of file TriggerMatcherToHLTDebug.cc.

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

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

Destructor.

Definition at line 194 of file TriggerMatcherToHLTDebug.cc.

194 {}

Member Function Documentation

void TriggerMatcherToHLTDebug::beginRun ( const edm::Run run,
const edm::EventSetup eventSetup 
)
overridevirtual

Reimplemented from edm::EDProducer.

Definition at line 405 of file TriggerMatcherToHLTDebug.cc.

References PropagateToMuon::init(), and l1matcher_.

405  {
406  l1matcher_.init(iSetup);
407 }
void init(const edm::EventSetup &iSetup)
Call this method at the beginning of each run, to initialize geometry, magnetic field and propagators...
void TriggerMatcherToHLTDebug::produce ( edm::Event event,
const edm::EventSetup eventSetup 
)
overridevirtual

Implements edm::EDProducer.

Definition at line 197 of file TriggerMatcherToHLTDebug.cc.

References SiPixelRawToDigiRegional_cfi::beamSpot, beamspotToken_, deltaR(), deltaR_, PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, reco::TrackBase::dxy(), reco::TrackBase::dz(), reco::TrackBase::error(), PV3DBase< T, PVType, FrameType >::eta(), reco::TrackBase::eta(), event(), PropagateToMuon::extrapolate(), TrajectoryStateOnSurface::globalPosition(), i, edm::HandleBase::id(), edm::Ref< C, T, F >::id(), 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(), HLT_25ns14e33_v1_cff::quality, seedMapToken_, edm::RefVector< C, T, F >::size(), tagToken_, theL2MuonsToken_, theL2SeedsToken_, theL3MuonsToken_, theL3SeedsToken_, and theL3TkTracksToken_.

197  {
198 
200  event.getByToken(tagToken_,muons);
201 
203  event.getByToken(l1Token_,L1Muons);
204 
206  event.getByToken(theL2SeedsToken_,L2Seeds);
207 
209  event.getByToken(theL2MuonsToken_,L2Muons);
210 
212  event.getByToken(theL3SeedsToken_,L3Seeds);
213 
215  event.getByToken(theL3TkTracksToken_,L3TkTracks);
216 
218  event.getByToken(theL3MuonsToken_,L3Muons);
219 
220  //beam spot
222  Handle<BeamSpot> recoBeamSpotHandle;
223  event.getByToken(beamspotToken_,recoBeamSpotHandle);
224  beamSpot = *recoBeamSpotHandle;
225 
226  //new for the MAP!!!!
227  edm::Handle<SeedMap> seedMapHandle;
228  event.getByToken(seedMapToken_, seedMapHandle);
229 
230 
231  size_t nmu = muons->size();
232  std::vector<int> propagatesToM2(nmu), hasL1Particle(nmu), hasL1Filtered(nmu);
233  std::vector<int> hasL2Seed(nmu), hasL2Muon(nmu), hasL2MuonFiltered(nmu);
234  std::vector<int> hasL3Seed(nmu), hasL3Track(nmu), hasL3TrackFiltered(nmu), hasL3Muon(nmu), hasL3MuonFiltered(nmu);
235  std::vector<reco::CandidatePtr> l1ptr(nmu), l2ptr(nmu), l3ptr(nmu);
236 
237  for (size_t i = 0; i < nmu; ++i) {
238  const reco::Muon &mu = (*muons)[i];
239 
240  // Propagate to muon station (using the L1 tool)
242  if (!stateAtMB2.isValid()) continue;
243  propagatesToM2[i] = 1;
244 
245  double etaTk = stateAtMB2.globalPosition().eta();
246  double phiTk = stateAtMB2.globalPosition().phi();
247  l1extra::L1MuonParticleCollection::const_iterator it;
248  vector<l1extra::L1MuonParticleRef>::const_iterator itMu3;
249  L2MuonTrajectorySeedCollection::const_iterator iSeed;
250  L3MuonTrajectorySeedCollection::const_iterator iSeedL3;
251  RecoChargedCandidateCollection::const_iterator iL2Muon;
252  reco::TrackCollection::const_iterator tktrackL3;
253  RecoChargedCandidateCollection::const_iterator iL3Muon;
254 
255  reco::CandidatePtr thisL1, thisL2, thisL3;
256  for(it = L1Muons->begin(); it != L1Muons->end(); ++it) {
257 
258  const L1MuGMTExtendedCand muonCand = (*it).gmtMuonCand();
259  unsigned int quality = muonCand.quality();
260 
261  double L1phi =(*it).phi();
262  double L1eta =(*it).eta();
263  double L1pt =(*it).pt();
264  double dR=deltaR(etaTk,phiTk,L1eta,L1phi);
265 
266  //CONDIZIONE-> CE NE E' UNA ASSOCIATA?
267  if (dR >= deltaR_) continue;
268  thisL1 = reco::CandidatePtr(L1Muons, it - L1Muons->begin());
269  if (!hasL1Particle[i]) l1ptr[i] = thisL1; // if nobody reached L1 before, then we're the best L1 found up to now.
270  hasL1Particle[i]++;
271 
272  if ((quality <= 3) || (L1pt<7)) continue;
273  if (!hasL1Filtered[i]) l1ptr[i] = thisL1; // if nobody reached L1 before, then we're the best L1 found up to now.
274  hasL1Filtered[i]++;
275 
276  if(!L2Seeds.isValid()) continue;
277  //LOOP SULLA COLLEZIONE DEI SEED
278  for( iSeed = L2Seeds->begin(); iSeed != L2Seeds->end(); ++iSeed) {
279 
280  l1extra::L1MuonParticleRef l1FromSeed = iSeed->l1Particle();
281  if (l1FromSeed.id() != L1Muons.id()) throw cms::Exception("CorruptData") << "You're using a different L1 collection than the one used by L2 seeds.\n";
282  if (l1FromSeed.key() != thisL1.key()) continue;
283  if (!hasL2Seed[i]) l1ptr[i] = thisL1; // if nobody reached here before, we're the best L1
284  hasL2Seed[i]++;
285 
286  if(!L2Muons.isValid()) continue;
287  //LOOP SULLA COLLEZIONE L2MUON
288  for( iL2Muon = L2Muons->begin(); iL2Muon != L2Muons->end(); ++iL2Muon) {
289 
290 
291  //MI FACCIO DARE REF E GUARDO SE E' UGUALE AL L2SEED ASSOCIATO
292  //BEFORE THE ASSOCIATION MAP!!!!!
293  //edm::Ref<L2MuonTrajectorySeedCollection> l2seedRef = iL2Muon->track()->seedRef().castTo<edm::Ref<L2MuonTrajectorySeedCollection> >();
294  //l1extra::L1MuonParticleRef l1FromL2 = l2seedRef->l1Particle();
295 
296  //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";
297  //if (l1FromL2 != l1FromSeed) continue;
298 
299  //AFTER THE ASSOCIATION MAP
300  const edm::RefVector<L2MuonTrajectorySeedCollection>& seeds = (*seedMapHandle)[iL2Muon->track()->seedRef().castTo<edm::Ref<L2MuonTrajectorySeedCollection> >()];
301  // bool isTriggered = false;
302  for(size_t jjj=0; jjj<seeds.size(); jjj++){
303 
304  if(seeds[jjj]->l1Particle()!= l1FromSeed) continue;
305 
306  }
307 
308 
309  thisL2 = reco::CandidatePtr(L2Muons, iL2Muon - L2Muons->begin()) ;
310  if (!hasL2Muon[i]) { l1ptr[i] = thisL1; l2ptr[i] = thisL2; } // if nobody reached here before, we're the best L1 and L2)
311  hasL2Muon[i]++;
312 
313  LogTrace(metname) <<"L2MUON TROVATO!"<<endl;
314  const reco::Track & L2Track = *iL2Muon->track();
315  double Eta_L2= L2Track.eta();
316  double Pt_L2= L2Track.pt();
317  int nValidHits_L2= L2Track.numberOfValidHits();
318  double BSPos_L2 = L2Track.dxy(beamSpot.position());
319  double dz_L2 =L2Track.dz();
320  double err0_L2 = L2Track.error(0);
321  double abspar0_L2 = fabs(L2Track.parameter(0));
322  double ptLx_L2 = Pt_L2;
323  if (abspar0_L2>0) ptLx_L2 += nsigma_Pt_L2*err0_L2/abspar0_L2*Pt_L2;
324 
325  //GUARDO SE L2MUON ASSOCIATO AVREBBE PASSATO IL FILTRO
326  bool passFilter = (((fabs(Eta_L2))<=max_Eta_L2)&&(nValidHits_L2>=min_Nhits_L2)&&((fabs(BSPos_L2))<=max_Dr_L2)&&((fabs(dz_L2))<=max_Dz_L2)&&(ptLx_L2>=min_Pt_L2));
327  if (!passFilter) continue;
328  if (!hasL2MuonFiltered[i]) { l1ptr[i] = thisL1; l2ptr[i] = thisL2; } // if nobody reached here before, we're the best L1 and L2)
329  hasL2MuonFiltered[i]++;
330 
331  const reco::TrackRef L2FilteredRef = iL2Muon->track();
332 
333  //########L3 PART##############
334  if (!L3Seeds.isValid()) continue;
335  for (iSeedL3 = L3Seeds->begin(); iSeedL3!= L3Seeds->end(); ++iSeedL3){
336 
337  TrackRef staTrack = iSeedL3->l2Track();
338  if (staTrack!=L2FilteredRef) continue;
339  if (!hasL3Seed[i]) { l1ptr[i] = thisL1; l2ptr[i] = thisL2; } // if nobody reached here before, we're the best L1 and L2)
340  hasL3Seed[i]++;
341 
342  if (!L3TkTracks.isValid()) continue;
343  for (tktrackL3 = L3TkTracks->begin(); tktrackL3!= L3TkTracks->end(); ++tktrackL3){
344 
346  TrackRef staTrack2 = l3seedRef->l2Track();
347 
348  if (staTrack2!=L2FilteredRef) continue;
349  if (!hasL3Track[i]) { l1ptr[i] = thisL1; l2ptr[i] = thisL2; } // if nobody reached here before, we're the best L1 and L2)
350  hasL3Track[i]++;
351 
352  if (!L3Muons.isValid()) continue;
353  for (iL3Muon = L3Muons->begin(); iL3Muon != L3Muons->end(); ++iL3Muon) {
354 
355  edm::Ref<L3MuonTrajectorySeedCollection> l3seedRef2 = iL3Muon->track()->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
356  TrackRef staTrack3 = l3seedRef2->l2Track();
357 
358  if (staTrack3!=L2FilteredRef) continue;
359  thisL3 = reco::CandidatePtr(L3Muons, iL3Muon - L3Muons->begin());
360 
361  if (!hasL3Muon[i]) { l1ptr[i] = thisL1; l2ptr[i] = thisL2; l3ptr[i] = thisL3; } // if nobody reached here before, we're the best L1, L2, L3
362  hasL3Muon[i]++;
363 
364  const reco::Track &L3Track = *iL3Muon->track();
365  double Eta_L3= L3Track.eta();
366  double Pt_L3= L3Track.pt();
367  int nValidHits_L3= L3Track.numberOfValidHits();
368  double BSPos_L3 = L3Track.dxy(beamSpot.position());
369  double dz_L3 =L3Track.dz();
370  double err0_L3 = L3Track.error(0);
371  double abspar0_L3 = fabs(L3Track.parameter(0));
372  double ptLx_L3 = Pt_L3;
373 
374  if (abspar0_L3>0) ptLx_L3 += nsigma_Pt_L3*err0_L3/abspar0_L3*Pt_L3;
375 
376  if(((fabs(Eta_L3))<=max_Eta_L3)&&(nValidHits_L3>=min_Nhits_L3)&&((fabs(BSPos_L3))<=max_Dr_L3)&&((fabs(dz_L3))<=max_Dz_L3)&&(ptLx_L3>=min_Pt_L3)){
377 
378  if (!hasL3MuonFiltered[i]) { l1ptr[i] = thisL1; l2ptr[i] = thisL2; l3ptr[i] = thisL3; } // if nobody reached here before, we're the best L1, L2, L3
379  hasL3MuonFiltered[i]++;
380 
381  }//L3MUON FILTERED ASSOCIATO TROVATO
382  }//L3MUON LOOP
383  }// L3 TRACKS
384  }// L3 SEEDS
385  }//T L2 MUONS
386  }// L2 SEEDS
387  }//L1 MUONS
388  } // RECO MUONS
389  storeValueMap<int>(event, muons, propagatesToM2, "propagatesToM2");
390  storeValueMap<int>(event, muons, hasL1Particle, "hasL1Particle");
391  storeValueMap<int>(event, muons, hasL1Filtered, "hasL1Filtered");
392  storeValueMap<int>(event, muons, hasL2Seed, "hasL2Seed");
393  storeValueMap<int>(event, muons, hasL2Muon, "hasL2Muon");
394  storeValueMap<int>(event, muons, hasL2MuonFiltered, "hasL2MuonFiltered");
395  storeValueMap<int>(event, muons, hasL3Seed, "hasL3Seed");
396  storeValueMap<int>(event, muons, hasL3Track, "hasL3Track");
397  storeValueMap<int>(event, muons, hasL3Muon, "hasL3Muon");
398  storeValueMap<int>(event, muons, hasL3MuonFiltered, "hasL3MuonFiltered");
399  storeValueMap<reco::CandidatePtr>(event, muons, l1ptr, "l1Candidate");
400  storeValueMap<reco::CandidatePtr>(event, muons, l2ptr, "l2Candidate");
401  storeValueMap<reco::CandidatePtr>(event, muons, l3ptr, "l3Candidate");
402 } // METHOD
int i
Definition: DBlmapReader.cc:9
TrajectoryStateOnSurface extrapolate(const reco::Track &tk) const
Extrapolate reco::Track to the muon station 2, return an invalid TSOS if it fails.
tuple L3Muons
Definition: L3Muons_cfi.py:30
key_type key() const
Definition: Ptr.h:186
ProductID id() const
Definition: HandleBase.cc:15
edm::EDGetTokenT< l1extra::L1MuonParticleCollection > l1Token_
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
edm::EDGetTokenT< L3MuonTrajectorySeedCollection > theL3SeedsToken_
edm::EDGetTokenT< L2MuonTrajectorySeedCollection > theL2SeedsToken_
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > theL2MuonsToken_
GlobalPoint globalPosition() const
key_type key() const
Accessor for product key.
Definition: Ref.h:264
tuple L2Muons
Definition: L2Muons_cfi.py:6
edm::EDGetTokenT< reco::RecoChargedCandidateCollection > theL3MuonsToken_
ProductID id() const
Accessor for product ID.
Definition: Ref.h:258
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:638
double pt() const
track transverse momentum
Definition: TrackBase.h:608
double error(int i) const
error on specified element
Definition: TrackBase.h:738
edm::EDGetTokenT< SeedMap > seedMapToken_
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:807
const int mu
Definition: Constants.h:22
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool isValid() const
Definition: HandleBase.h:75
#define LogTrace(id)
double parameter(int i) const
i-th parameter ( i = 0, ... 4 )
Definition: TrackBase.h:726
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:596
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
edm::EDGetTokenT< edm::View< reco::Muon > > tagToken_
unsigned int quality() const
get quality
Definition: L1MuGMTCand.h:93
edm::Ptr< Candidate > CandidatePtr
persistent reference to an object in a collection of Candidate objects
Definition: CandidateFwd.h:25
T eta() const
Definition: PV3DBase.h:76
tuple muons
Definition: patZpeak.py:38
edm::EDGetTokenT< reco::BeamSpot > beamspotToken_
edm::EDGetTokenT< reco::TrackCollection > theL3TkTracksToken_
size_type size() const
Size of the RefVector.
Definition: RefVector.h:99
const Point & position() const
position
Definition: BeamSpot.h:62
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:578
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 411 of file TriggerMatcherToHLTDebug.cc.

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

414  {
415  using namespace edm; using namespace std;
416  auto_ptr<ValueMap<T> > valMap(new ValueMap<T>());
417  typename edm::ValueMap<T>::Filler filler(*valMap);
418  filler.insert(handle, values.begin(), values.end());
419  filler.fill();
420  iEvent.put(valMap, label);
421 }
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120

Member Data Documentation

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

Definition at line 101 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::deltaR_
private

Definition at line 95 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

PropagateToMuon TriggerMatcherToHLTDebug::l1matcher_
private

Definition at line 90 of file TriggerMatcherToHLTDebug.cc.

Referenced by beginRun(), and produce().

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

Definition at line 89 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Dr_L2
private

Definition at line 105 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Dr_L3
private

Definition at line 114 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Dz_L2
private

Definition at line 106 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Dz_L3
private

Definition at line 115 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Eta_L2
private

Definition at line 103 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Eta_L3
private

Definition at line 112 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

std::string TriggerMatcherToHLTDebug::metname
private

Definition at line 92 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

int TriggerMatcherToHLTDebug::min_N_L2
private

Definition at line 102 of file TriggerMatcherToHLTDebug.cc.

int TriggerMatcherToHLTDebug::min_N_L3
private

Definition at line 111 of file TriggerMatcherToHLTDebug.cc.

int TriggerMatcherToHLTDebug::min_Nhits_L2
private

Definition at line 104 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

int TriggerMatcherToHLTDebug::min_Nhits_L3
private

Definition at line 113 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::min_Pt_L2
private

Definition at line 107 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::min_Pt_L3
private

Definition at line 116 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

int TriggerMatcherToHLTDebug::minL1Quality_
private

Definition at line 98 of file TriggerMatcherToHLTDebug.cc.

double TriggerMatcherToHLTDebug::nsigma_Pt_L2
private

Definition at line 108 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::nsigma_Pt_L3
private

Definition at line 117 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

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

Definition at line 124 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

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

Definition at line 88 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

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

Definition at line 120 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

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

Definition at line 119 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

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

Definition at line 123 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

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

Definition at line 121 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

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

Definition at line 122 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().