CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MuonAlignmentPreFilter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: MuonAlignmentPreFilter
4 // Class: MuonAlignmentPreFilter
5 //
13 #include <memory>
14 
22 
29 
30 
32 {
33 public:
36 
37 private:
38  virtual void beginJob() override {}
39  virtual bool filter(edm::Event&, const edm::EventSetup&) override;
40  virtual void endJob() override {}
41 
42  // ----------member data ---------------------------
43 
45  double m_minTrackPt;
46  double m_minTrackP;
51  double m_minTrackEta;
52  double m_maxTrackEta;
53 };
54 
55 
57  : m_tracksTag(cfg.getParameter<edm::InputTag>("tracksTag"))
58  , m_minTrackPt(cfg.getParameter<double>("minTrackPt"))
59  , m_minTrackP(cfg.getParameter<double>("minTrackP"))
60  , m_allowTIDTEC(cfg.getParameter<bool>("allowTIDTEC"))
61  , m_minTrackerHits(cfg.getParameter<int>("minTrackerHits"))
62  , m_minDTHits(cfg.getParameter<int>("minDTHits"))
63  , m_minCSCHits(cfg.getParameter<int>("minCSCHits"))
64  , m_minTrackEta(cfg.getParameter<double>("minTrackEta"))
65  , m_maxTrackEta(cfg.getParameter<double>("maxTrackEta"))
66 {}
67 
68 
69 bool
71 {
73  iEvent.getByLabel(m_tracksTag, trackColl);
74 
75  // check if there's at least one interesting track:
76 
77  for (reco::TrackCollection::const_iterator it = trackColl->begin(); it != trackColl->end(); it++)
78  {
79  int tracker_numHits = 0;
80  bool contains_TIDTEC = false;
81  int dt_numHits = 0;
82  int csc_numHits = 0;
83 
84  const reco::Track* track = &(*it);
85 
86  if (track->pt() < m_minTrackPt || track->p() < m_minTrackP) continue;
87  if (track->eta() < m_minTrackEta || track->eta() > m_maxTrackEta ) continue;
88 
89  for (trackingRecHit_iterator hit = track->recHitsBegin(); hit != track->recHitsEnd(); ++hit)
90  {
91  DetId id = (*hit)->geographicalId();
92  if (id.det() == DetId::Tracker)
93  {
94  tracker_numHits++;
95  if (id.subdetId() == StripSubdetector::TID || id.subdetId() == StripSubdetector::TEC) contains_TIDTEC = true;
96  }
97 
98  if ( id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::DT ) dt_numHits++;
99  if ( id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::CSC ) csc_numHits++;
100  }
101 
102  if (( m_allowTIDTEC || !contains_TIDTEC ) &&
103  m_minTrackerHits <= tracker_numHits &&
104  ( m_minDTHits <= dt_numHits || m_minCSCHits <= csc_numHits ) ) return true;
105  }
106  return false;
107 }
108 
109 
110 //define this as a plug-in
double p() const
momentum vector magnitude
Definition: TrackBase.h:602
tuple cfg
Definition: looper.py:293
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void endJob() override
int iEvent
Definition: GenABIO.cc:230
static const int CSC
Definition: MuonSubdetId.h:13
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:638
virtual void beginJob() override
double pt() const
track transverse momentum
Definition: TrackBase.h:608
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
virtual bool filter(edm::Event &, const edm::EventSetup &) override
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
Definition: DetId.h:18
static const int DT
Definition: MuonSubdetId.h:12
MuonAlignmentPreFilter(const edm::ParameterSet &)
TrackingRecHitCollection::base::const_iterator trackingRecHit_iterator
iterator over a vector of reference to TrackingRecHit in the same collection
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:109