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 ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 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
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) 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::InputTag beamspotTag_
 
double deltaR_
 
PropagateToMuon l1matcher_
 
edm::InputTag l1Tag_
 
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::InputTag seedMapTag_
 
edm::InputTag tagTag_
 
edm::InputTag theL2MuonsLabel
 
edm::InputTag theL2SeedsLabel
 
edm::InputTag theL3MuonsLabel
 
edm::InputTag theL3SeedsLabel
 
edm::InputTag theL3TkTracksLabel
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- 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::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- 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 74 of file TriggerMatcherToHLTDebug.cc.

Member Typedef Documentation

Definition at line 126 of file TriggerMatcherToHLTDebug.cc.

Constructor & Destructor Documentation

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

Definition at line 144 of file TriggerMatcherToHLTDebug.cc.

References edm::ParameterSet::getParameter(), metname, theL2MuonsLabel, theL2SeedsLabel, theL3MuonsLabel, theL3SeedsLabel, and theL3TkTracksLabel.

144  :
145  tagTag_(pset.getParameter<edm::InputTag>("tags")),
146  l1Tag_(pset.getParameter<edm::InputTag>("l1s")),
147  l1matcher_(pset.getParameter<edm::ParameterSet>("l1matcherConfig")),
148  deltaR_(pset.getParameter<double>("deltaR")),
149  minL1Quality_(pset.getParameter<int32_t>("MinL1Quality")),
150  beamspotTag_(pset.getParameter<edm::InputTag>("BeamSpotTag")),
151  min_N_L2(pset.getParameter<int> ("MinN_L2")),
152  max_Eta_L2(pset.getParameter<double> ("MaxEta_L2")),
153  min_Nhits_L2(pset.getParameter<int> ("MinNhits_L2")),
154  max_Dr_L2(pset.getParameter<double> ("MaxDr_L2")),
155  max_Dz_L2(pset.getParameter<double> ("MaxDz_L2")),
156  min_Pt_L2(pset.getParameter<double> ("MinPt_L2")),
157  nsigma_Pt_L2(pset.getParameter<double> ("NSigmaPt_L2")),
158  min_N_L3(pset.getParameter<int> ("MinN_L3")),
159  max_Eta_L3(pset.getParameter<double> ("MaxEta_L3")),
160  min_Nhits_L3(pset.getParameter<int> ("MinNhits_L3")),
161  max_Dr_L3(pset.getParameter<double> ("MaxDr_L3")),
162  max_Dz_L3(pset.getParameter<double> ("MaxDz_L3")),
163  min_Pt_L3(pset.getParameter<double> ("MinPt_L3")),
164  nsigma_Pt_L3(pset.getParameter<double> ("NSigmaPt_L3")),
165  seedMapTag_(pset.getParameter<edm::InputTag >("SeedMapTag"))
166 {
167 
168 
169  theL2SeedsLabel = pset.getParameter<InputTag>("L2Seeds_Collection");
170  theL2MuonsLabel = pset.getParameter<InputTag>("L2Muons_Collection");
171  theL3SeedsLabel = pset.getParameter<InputTag>("L3Seeds_Collection");
172  theL3TkTracksLabel = pset.getParameter<InputTag>("L3TkTracks_Collection");
173  theL3MuonsLabel = pset.getParameter<InputTag>("L3Muons_Collection");
174 
175 
176  metname = "TriggerMatcherToHLTDebug";
177 
178  produces<edm::ValueMap<int> > ("propagatesToM2");
179  produces<edm::ValueMap<int> > ("hasL1Particle");
180  produces<edm::ValueMap<int> > ("hasL1Filtered");
181  produces<edm::ValueMap<int> > ("hasL2Seed");
182  produces<edm::ValueMap<int> > ("hasL2Muon");
183  produces<edm::ValueMap<int> > ("hasL2MuonFiltered");
184  produces<edm::ValueMap<int> > ("hasL3Seed");
185  produces<edm::ValueMap<int> > ("hasL3Track");
186  produces<edm::ValueMap<int> > ("hasL3Muon");
187  produces<edm::ValueMap<int> > ("hasL3MuonFiltered");
188 
189  produces<edm::ValueMap<reco::CandidatePtr> > ("l1Candidate");
190  produces<edm::ValueMap<reco::CandidatePtr> > ("l2Candidate");
191  produces<edm::ValueMap<reco::CandidatePtr> > ("l3Candidate");
192 }
T getParameter(std::string const &) const
TriggerMatcherToHLTDebug::~TriggerMatcherToHLTDebug ( )
virtual

Destructor.

Definition at line 197 of file TriggerMatcherToHLTDebug.cc.

197 {}

Member Function Documentation

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

Reimplemented from edm::EDProducer.

Definition at line 408 of file TriggerMatcherToHLTDebug.cc.

References PropagateToMuon::init(), and l1matcher_.

408  {
409  l1matcher_.init(iSetup);
410 }
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 200 of file TriggerMatcherToHLTDebug.cc.

References SiPixelRawToDigiRegional_cfi::beamSpot, beamspotTag_, 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(), edm::HandleBase::isValid(), TrajectoryStateOnSurface::isValid(), edm::Ptr< T >::key(), edm::Ref< C, T, F >::key(), l1matcher_, l1Tag_, 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(), seedMapTag_, edm::RefVector< C, T, F >::size(), tagTag_, theL2MuonsLabel, theL2SeedsLabel, theL3MuonsLabel, theL3SeedsLabel, and theL3TkTracksLabel.

200  {
201 
203  event.getByLabel(tagTag_,muons);
204 
206  event.getByLabel(l1Tag_,L1Muons);
207 
209  event.getByLabel(theL2SeedsLabel,L2Seeds);
210 
212  event.getByLabel(theL2MuonsLabel,L2Muons);
213 
215  event.getByLabel(theL3SeedsLabel,L3Seeds);
216 
218  event.getByLabel(theL3TkTracksLabel,L3TkTracks);
219 
221  event.getByLabel(theL3MuonsLabel,L3Muons);
222 
223  //beam spot
225  Handle<BeamSpot> recoBeamSpotHandle;
226  event.getByLabel(beamspotTag_,recoBeamSpotHandle);
227  beamSpot = *recoBeamSpotHandle;
228 
229  //new for the MAP!!!!
230  edm::Handle<SeedMap> seedMapHandle;
231  event.getByLabel(seedMapTag_, seedMapHandle);
232 
233 
234  size_t nmu = muons->size();
235  std::vector<int> propagatesToM2(nmu), hasL1Particle(nmu), hasL1Filtered(nmu);
236  std::vector<int> hasL2Seed(nmu), hasL2Muon(nmu), hasL2MuonFiltered(nmu);
237  std::vector<int> hasL3Seed(nmu), hasL3Track(nmu), hasL3TrackFiltered(nmu), hasL3Muon(nmu), hasL3MuonFiltered(nmu);
238  std::vector<reco::CandidatePtr> l1ptr(nmu), l2ptr(nmu), l3ptr(nmu);
239 
240  for (size_t i = 0; i < nmu; ++i) {
241  const reco::Muon &mu = (*muons)[i];
242 
243  // Propagate to muon station (using the L1 tool)
245  if (!stateAtMB2.isValid()) continue;
246  propagatesToM2[i] = 1;
247 
248  double etaTk = stateAtMB2.globalPosition().eta();
249  double phiTk = stateAtMB2.globalPosition().phi();
250  l1extra::L1MuonParticleCollection::const_iterator it;
251  vector<l1extra::L1MuonParticleRef>::const_iterator itMu3;
252  L2MuonTrajectorySeedCollection::const_iterator iSeed;
253  L3MuonTrajectorySeedCollection::const_iterator iSeedL3;
254  RecoChargedCandidateCollection::const_iterator iL2Muon;
255  reco::TrackCollection::const_iterator tktrackL3;
256  RecoChargedCandidateCollection::const_iterator iL3Muon;
257 
258  reco::CandidatePtr thisL1, thisL2, thisL3;
259  for(it = L1Muons->begin(); it != L1Muons->end(); ++it) {
260 
261  const L1MuGMTExtendedCand muonCand = (*it).gmtMuonCand();
262  unsigned int quality = muonCand.quality();
263 
264  double L1phi =(*it).phi();
265  double L1eta =(*it).eta();
266  double L1pt =(*it).pt();
267  double dR=deltaR(etaTk,phiTk,L1eta,L1phi);
268 
269  //CONDIZIONE-> CE NE E' UNA ASSOCIATA?
270  if (dR >= deltaR_) continue;
271  thisL1 = reco::CandidatePtr(L1Muons, it - L1Muons->begin());
272  if (!hasL1Particle[i]) l1ptr[i] = thisL1; // if nobody reached L1 before, then we're the best L1 found up to now.
273  hasL1Particle[i]++;
274 
275  if ((quality <= 3) || (L1pt<7)) continue;
276  if (!hasL1Filtered[i]) l1ptr[i] = thisL1; // if nobody reached L1 before, then we're the best L1 found up to now.
277  hasL1Filtered[i]++;
278 
279  if(!L2Seeds.isValid()) continue;
280  //LOOP SULLA COLLEZIONE DEI SEED
281  for( iSeed = L2Seeds->begin(); iSeed != L2Seeds->end(); ++iSeed) {
282 
283  l1extra::L1MuonParticleRef l1FromSeed = iSeed->l1Particle();
284  if (l1FromSeed.id() != L1Muons.id()) throw cms::Exception("CorruptData") << "You're using a different L1 collection than the one used by L2 seeds.\n";
285  if (l1FromSeed.key() != thisL1.key()) continue;
286  if (!hasL2Seed[i]) l1ptr[i] = thisL1; // if nobody reached here before, we're the best L1
287  hasL2Seed[i]++;
288 
289  if(!L2Muons.isValid()) continue;
290  //LOOP SULLA COLLEZIONE L2MUON
291  for( iL2Muon = L2Muons->begin(); iL2Muon != L2Muons->end(); ++iL2Muon) {
292 
293 
294  //MI FACCIO DARE REF E GUARDO SE E' UGUALE AL L2SEED ASSOCIATO
295  //BEFORE THE ASSOCIATION MAP!!!!!
296  //edm::Ref<L2MuonTrajectorySeedCollection> l2seedRef = iL2Muon->track()->seedRef().castTo<edm::Ref<L2MuonTrajectorySeedCollection> >();
297  //l1extra::L1MuonParticleRef l1FromL2 = l2seedRef->l1Particle();
298 
299  //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";
300  //if (l1FromL2 != l1FromSeed) continue;
301 
302  //AFTER THE ASSOCIATION MAP
303  const edm::RefVector<L2MuonTrajectorySeedCollection>& seeds = (*seedMapHandle)[iL2Muon->track()->seedRef().castTo<edm::Ref<L2MuonTrajectorySeedCollection> >()];
304  // bool isTriggered = false;
305  for(size_t jjj=0; jjj<seeds.size(); jjj++){
306 
307  if(seeds[jjj]->l1Particle()!= l1FromSeed) continue;
308 
309  }
310 
311 
312  thisL2 = reco::CandidatePtr(L2Muons, iL2Muon - L2Muons->begin()) ;
313  if (!hasL2Muon[i]) { l1ptr[i] = thisL1; l2ptr[i] = thisL2; } // if nobody reached here before, we're the best L1 and L2)
314  hasL2Muon[i]++;
315 
316  LogTrace(metname) <<"L2MUON TROVATO!"<<endl;
317  const reco::Track & L2Track = *iL2Muon->track();
318  double Eta_L2= L2Track.eta();
319  double Pt_L2= L2Track.pt();
320  int nValidHits_L2= L2Track.numberOfValidHits();
321  double BSPos_L2 = L2Track.dxy(beamSpot.position());
322  double dz_L2 =L2Track.dz();
323  double err0_L2 = L2Track.error(0);
324  double abspar0_L2 = fabs(L2Track.parameter(0));
325  double ptLx_L2 = Pt_L2;
326  if (abspar0_L2>0) ptLx_L2 += nsigma_Pt_L2*err0_L2/abspar0_L2*Pt_L2;
327 
328  //GUARDO SE L2MUON ASSOCIATO AVREBBE PASSATO IL FILTRO
329  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));
330  if (!passFilter) continue;
331  if (!hasL2MuonFiltered[i]) { l1ptr[i] = thisL1; l2ptr[i] = thisL2; } // if nobody reached here before, we're the best L1 and L2)
332  hasL2MuonFiltered[i]++;
333 
334  const reco::TrackRef L2FilteredRef = iL2Muon->track();
335 
336  //########L3 PART##############
337  if (!L3Seeds.isValid()) continue;
338  for (iSeedL3 = L3Seeds->begin(); iSeedL3!= L3Seeds->end(); ++iSeedL3){
339 
340  TrackRef staTrack = iSeedL3->l2Track();
341  if (staTrack!=L2FilteredRef) continue;
342  if (!hasL3Seed[i]) { l1ptr[i] = thisL1; l2ptr[i] = thisL2; } // if nobody reached here before, we're the best L1 and L2)
343  hasL3Seed[i]++;
344 
345  if (!L3TkTracks.isValid()) continue;
346  for (tktrackL3 = L3TkTracks->begin(); tktrackL3!= L3TkTracks->end(); ++tktrackL3){
347 
349  TrackRef staTrack2 = l3seedRef->l2Track();
350 
351  if (staTrack2!=L2FilteredRef) continue;
352  if (!hasL3Track[i]) { l1ptr[i] = thisL1; l2ptr[i] = thisL2; } // if nobody reached here before, we're the best L1 and L2)
353  hasL3Track[i]++;
354 
355  if (!L3Muons.isValid()) continue;
356  for (iL3Muon = L3Muons->begin(); iL3Muon != L3Muons->end(); ++iL3Muon) {
357 
358  edm::Ref<L3MuonTrajectorySeedCollection> l3seedRef2 = iL3Muon->track()->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
359  TrackRef staTrack3 = l3seedRef2->l2Track();
360 
361  if (staTrack3!=L2FilteredRef) continue;
362  thisL3 = reco::CandidatePtr(L3Muons, iL3Muon - L3Muons->begin());
363 
364  if (!hasL3Muon[i]) { l1ptr[i] = thisL1; l2ptr[i] = thisL2; l3ptr[i] = thisL3; } // if nobody reached here before, we're the best L1, L2, L3
365  hasL3Muon[i]++;
366 
367  const reco::Track &L3Track = *iL3Muon->track();
368  double Eta_L3= L3Track.eta();
369  double Pt_L3= L3Track.pt();
370  int nValidHits_L3= L3Track.numberOfValidHits();
371  double BSPos_L3 = L3Track.dxy(beamSpot.position());
372  double dz_L3 =L3Track.dz();
373  double err0_L3 = L3Track.error(0);
374  double abspar0_L3 = fabs(L3Track.parameter(0));
375  double ptLx_L3 = Pt_L3;
376 
377  if (abspar0_L3>0) ptLx_L3 += nsigma_Pt_L3*err0_L3/abspar0_L3*Pt_L3;
378 
379  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)){
380 
381  if (!hasL3MuonFiltered[i]) { l1ptr[i] = thisL1; l2ptr[i] = thisL2; l3ptr[i] = thisL3; } // if nobody reached here before, we're the best L1, L2, L3
382  hasL3MuonFiltered[i]++;
383 
384  }//L3MUON FILTERED ASSOCIATO TROVATO
385  }//L3MUON LOOP
386  }// L3 TRACKS
387  }// L3 SEEDS
388  }//T L2 MUONS
389  }// L2 SEEDS
390  }//L1 MUONS
391  } // RECO MUONS
392  storeValueMap<int>(event, muons, propagatesToM2, "propagatesToM2");
393  storeValueMap<int>(event, muons, hasL1Particle, "hasL1Particle");
394  storeValueMap<int>(event, muons, hasL1Filtered, "hasL1Filtered");
395  storeValueMap<int>(event, muons, hasL2Seed, "hasL2Seed");
396  storeValueMap<int>(event, muons, hasL2Muon, "hasL2Muon");
397  storeValueMap<int>(event, muons, hasL2MuonFiltered, "hasL2MuonFiltered");
398  storeValueMap<int>(event, muons, hasL3Seed, "hasL3Seed");
399  storeValueMap<int>(event, muons, hasL3Track, "hasL3Track");
400  storeValueMap<int>(event, muons, hasL3Muon, "hasL3Muon");
401  storeValueMap<int>(event, muons, hasL3MuonFiltered, "hasL3MuonFiltered");
402  storeValueMap<reco::CandidatePtr>(event, muons, l1ptr, "l1Candidate");
403  storeValueMap<reco::CandidatePtr>(event, muons, l2ptr, "l2Candidate");
404  storeValueMap<reco::CandidatePtr>(event, muons, l3ptr, "l3Candidate");
405 } // 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
ProductID id() const
Definition: HandleBase.cc:15
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
GlobalPoint globalPosition() const
tuple L2Muons
Definition: L2Muons_cfi.py:6
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:141
double pt() const
track transverse momentum
Definition: TrackBase.h:131
double error(int i) const
error on specified element
Definition: TrackBase.h:189
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:232
const int mu
Definition: Constants.h:23
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:76
#define LogTrace(id)
double parameter(int i) const
i-th parameter ( i = 0, ... 4 )
Definition: TrackBase.h:185
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:127
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
unsigned int quality() const
get quality
Definition: L1MuGMTCand.h:95
key_type key() const
Definition: Ptr.h:169
edm::Ptr< Candidate > CandidatePtr
persistent reference to an object in a collection of Candidate objects
Definition: CandidateFwd.h:25
key_type key() const
Accessor for product key.
Definition: Ref.h:266
T eta() const
Definition: PV3DBase.h:76
tuple muons
Definition: patZpeak.py:38
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
const Point & position() const
position
Definition: BeamSpot.h:63
ProductID id() const
Accessor for product ID.
Definition: Ref.h:256
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:121
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 414 of file TriggerMatcherToHLTDebug.cc.

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

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

Member Data Documentation

edm::InputTag TriggerMatcherToHLTDebug::beamspotTag_
private

Definition at line 100 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::deltaR_
private

Definition at line 94 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

PropagateToMuon TriggerMatcherToHLTDebug::l1matcher_
private

Definition at line 89 of file TriggerMatcherToHLTDebug.cc.

Referenced by beginRun(), and produce().

edm::InputTag TriggerMatcherToHLTDebug::l1Tag_
private

Definition at line 88 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Dr_L2
private

Definition at line 104 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Dr_L3
private

Definition at line 113 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Dz_L2
private

Definition at line 105 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Dz_L3
private

Definition at line 114 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Eta_L2
private

Definition at line 102 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::max_Eta_L3
private

Definition at line 111 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

std::string TriggerMatcherToHLTDebug::metname
private

Definition at line 91 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

int TriggerMatcherToHLTDebug::min_N_L2
private

Definition at line 101 of file TriggerMatcherToHLTDebug.cc.

int TriggerMatcherToHLTDebug::min_N_L3
private

Definition at line 110 of file TriggerMatcherToHLTDebug.cc.

int TriggerMatcherToHLTDebug::min_Nhits_L2
private

Definition at line 103 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

int TriggerMatcherToHLTDebug::min_Nhits_L3
private

Definition at line 112 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::min_Pt_L2
private

Definition at line 106 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::min_Pt_L3
private

Definition at line 115 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

int TriggerMatcherToHLTDebug::minL1Quality_
private

Definition at line 97 of file TriggerMatcherToHLTDebug.cc.

double TriggerMatcherToHLTDebug::nsigma_Pt_L2
private

Definition at line 107 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

double TriggerMatcherToHLTDebug::nsigma_Pt_L3
private

Definition at line 116 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

edm::InputTag TriggerMatcherToHLTDebug::seedMapTag_
private

Definition at line 123 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

edm::InputTag TriggerMatcherToHLTDebug::tagTag_
private

Definition at line 88 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce().

edm::InputTag TriggerMatcherToHLTDebug::theL2MuonsLabel
private

Definition at line 119 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

edm::InputTag TriggerMatcherToHLTDebug::theL2SeedsLabel
private

Definition at line 118 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

edm::InputTag TriggerMatcherToHLTDebug::theL3MuonsLabel
private

Definition at line 122 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

edm::InputTag TriggerMatcherToHLTDebug::theL3SeedsLabel
private

Definition at line 120 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().

edm::InputTag TriggerMatcherToHLTDebug::theL3TkTracksLabel
private

Definition at line 121 of file TriggerMatcherToHLTDebug.cc.

Referenced by produce(), and TriggerMatcherToHLTDebug().