CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
DTSegmentSelector.cc
Go to the documentation of this file.
2 
8 
15 
17  : muonTags_(pset.getParameter<edm::InputTag>("Muons")),
18  checkNoisyChannels_(pset.getParameter<bool>("checkNoisyChannels")),
19  minHitsPhi_(pset.getParameter<int>("minHitsPhi")),
20  minHitsZ_(pset.getParameter<int>("minHitsZ")),
21  maxChi2_(pset.getParameter<double>("maxChi2")),
22  maxAnglePhi_(pset.getParameter<double>("maxAnglePhi")),
23  maxAngleZ_(pset.getParameter<double>("maxAngleZ")) {
26 }
27 
29  edm::Event const& event,
30  edm::EventSetup const& setup) {
31  bool result = true;
32 
33  /* get the muon collection if one is specified
34  (not specifying a muon collection switches off muon matching */
35  if (!muonTags_.label().empty()) {
37  event.getByToken(muonToken_, muonH);
38  const std::vector<reco::Muon>* muons = muonH.product();
39  //std::cout << " Muon collection size: " << muons->size() << std::endl;
40  if (muons->empty())
41  return false;
42 
43  DTChamberId segId(segment.rawId());
44 
45  bool matched = false;
46  for (auto& imuon : *muons)
47  for (const auto& ch : imuon.matches()) {
48  DetId chId(ch.id.rawId());
49  if (chId != segId)
50  continue;
51  if (imuon.pt() < 15 || !imuon.isGlobalMuon())
52  continue;
53 
54  int nsegs = ch.segmentMatches.size();
55  if (!nsegs)
56  continue;
57 
58  LocalPoint posHit = segment.localPosition();
59  float dx = (posHit.x() ? posHit.x() - ch.x : 0);
60  float dy = (posHit.y() ? posHit.y() - ch.y : 0);
61  float dr = sqrt(dx * dx + dy * dy);
62  if (dr < 5)
63  matched = true;
64  }
65 
66  if (!matched)
67  result = false;
68  }
69 
72  statusMap = setup.getHandle(theStatusMapToken_);
73 
74  // Get the Phi 2D segment
75  int nPhiHits = -1;
76  bool segmentNoisyPhi = false;
77  if (segment.hasPhi()) {
78  const DTChamberRecSegment2D* phiSeg = segment.phiSegment(); // phiSeg lives in the chamber RF
79  std::vector<DTRecHit1D> phiRecHits = phiSeg->specificRecHits();
80  nPhiHits = phiRecHits.size();
82  segmentNoisyPhi = checkNoisySegment(statusMap, phiRecHits);
83  //} else result = false;
84  }
85  // Get the Theta 2D segment
86  int nZHits = -1;
87  bool segmentNoisyZ = false;
88  if (segment.hasZed()) {
89  const DTSLRecSegment2D* zSeg = segment.zSegment(); // zSeg lives in the SL RF
90  std::vector<DTRecHit1D> zRecHits = zSeg->specificRecHits();
91  nZHits = zRecHits.size();
93  segmentNoisyZ = checkNoisySegment(statusMap, zRecHits);
94  }
95 
96  // Segment selection
97  // Discard segment if it has a noisy cell
98  if (segmentNoisyPhi || segmentNoisyZ)
99  result = false;
100 
101  // 2D-segment number of hits
102  if (segment.hasPhi() && nPhiHits < minHitsPhi_)
103  result = false;
104 
105  if (segment.hasZed() && nZHits < minHitsZ_)
106  result = false;
107 
108  // Segment chi2
109  double chiSquare = segment.chi2() / segment.degreesOfFreedom();
110  if (chiSquare > maxChi2_)
111  result = false;
112 
113  // Segment angle
114  LocalVector segment4DLocalDir = segment.localDirection();
115  double angleZ = fabs(atan(segment4DLocalDir.y() / segment4DLocalDir.z()) * 180. / Geom::pi());
116  if (angleZ > maxAngleZ_)
117  result = false;
118 
119  double anglePhi = fabs(atan(segment4DLocalDir.x() / segment4DLocalDir.z()) * 180. / Geom::pi());
120  if (anglePhi > maxAnglePhi_)
121  result = false;
122 
123  return result;
124 }
125 
127  std::vector<DTRecHit1D> const& dtHits) {
128  bool segmentNoisy = false;
129 
130  std::vector<DTRecHit1D>::const_iterator dtHit = dtHits.begin();
131  std::vector<DTRecHit1D>::const_iterator dtHits_end = dtHits.end();
132  for (; dtHit != dtHits_end; ++dtHit) {
133  DTWireId wireId = dtHit->wireId();
134  // Check for noisy channels to skip them
135  bool isNoisy = false, isFEMasked = false, isTDCMasked = false, isTrigMask = false, isDead = false, isNohv = false;
136  statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
137  if (isNoisy) {
138  LogTrace("Calibration") << "Wire: " << wireId << " is noisy, skipping!";
139  segmentNoisy = true;
140  break;
141  }
142  }
143  return segmentNoisy;
144 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
DTSegmentSelector(edm::ParameterSet const &pset, edm::ConsumesCollector &iC)
const DTChamberRecSegment2D * phiSegment() const
The superPhi segment: 0 if no phi projection available.
LocalVector localDirection() const override
Local direction in Chamber frame.
T y() const
Definition: PV3DBase.h:60
LocalPoint localPosition() const override
Local position in Chamber frame.
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
#define LogTrace(id)
tuple result
Definition: mps_fire.py:311
bool operator()(DTRecSegment4D const &, edm::Event const &, edm::EventSetup const &)
T sqrt(T t)
Definition: SSEVec.h:19
T z() const
Definition: PV3DBase.h:61
edm::ESGetToken< DTStatusFlag, DTStatusFlagRcd > theStatusMapToken_
bool hasPhi() const
Does it have the Phi projection?
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
bool hasZed() const
Does it have the Z projection?
Definition: DetId.h:17
const DTSLRecSegment2D * zSegment() const
The Z segment: 0 if not zed projection available.
edm::EDGetTokenT< reco::MuonCollection > muonToken_
T const * product() const
Definition: Handle.h:70
int degreesOfFreedom() const override
Degrees of freedom of the segment fit.
bool checkNoisySegment(edm::ESHandle< DTStatusFlag > const &, std::vector< DTRecHit1D > const &)
std::string const & label() const
Definition: InputTag.h:36
tuple muons
Definition: patZpeak.py:41
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:151
constexpr double pi()
Definition: Pi.h:31
T x() const
Definition: PV3DBase.h:59
edm::InputTag muonTags_
id_type rawId() const
double chi2() const override
Chi2 of the segment fit.