CMS 3D CMS Logo

MuonSegmentMatcher.cc
Go to the documentation of this file.
1 // -*- C++ -*- // // Package: MuonSegmentMatcher // Class: MuonSegmentMatcher // /**\class MuonSegmentMatcher MuonSegmentMatcher.cc
2 // Description: <one line class summary>
3 // Implementation:
4 // <Notes on implementation>
5 //*/
6 //
7 // Original Author: Alan Tua
8 // Created: Wed Jul 9 21:40:17 CEST 2008
9 //
10 //
11 
13 
14 // user include files
17 
19 
22 
25 
31 
32 // system include files
33 #include <vector>
34 #include <iostream>
35 
36 class MuonSegmentMatcher;
37 
38 using namespace std;
39 
40 // constructors and destructor
41 
43  : DTSegmentTags_(matchParameters.getParameter<edm::InputTag>("DTsegments")),
44  CSCSegmentTags_(matchParameters.getParameter<edm::InputTag>("CSCsegments")),
45  dtRadius_(matchParameters.getParameter<double>("DTradius")),
46  dtTightMatch(matchParameters.getParameter<bool>("TightMatchDT")),
47  cscTightMatch(matchParameters.getParameter<bool>("TightMatchCSC")) {
48  if (matchParameters.existsAs<edm::InputTag>("RPChits")) {
49  RPCHitTags_ = matchParameters.getParameter<edm::InputTag>("RPChits");
50  } else
51  RPCHitTags_ = edm::InputTag("hltRpcRecHits");
52 
56 }
57 
59 
60 // member functions
61 vector<const DTRecSegment4D*> MuonSegmentMatcher::matchDT(const reco::Track& muon, const edm::Event& event) {
62  using namespace edm;
63 
65  event.getByToken(dtRecHitsToken, dtRecHits);
66 
67  vector<const DTRecSegment4D*> pointerTo4DSegments;
68 
69  std::vector<TrackingRecHit const*> dtHits;
70 
71  bool segments = false;
72 
73  // Loop and select DT recHits
74  for (auto const& hit : muon.recHits()) {
75  if (!hit->isValid())
76  continue;
77  if (hit->geographicalId().det() != DetId::Muon)
78  continue;
79  if (hit->geographicalId().subdetId() != MuonSubdetId::DT)
80  continue;
81  if (!hit->recHits().empty())
82  if ((*hit->recHits().begin())->recHits().size() > 1)
83  segments = true;
84  dtHits.push_back(hit);
85  }
86 
87  // cout << "Muon DT hits found: " << dtHits.size() << " segments " << segments << endl;
88 
89  double PhiCutParameter = dtRadius_;
90  double ZCutParameter = dtRadius_;
91  double matchRatioZ = 0;
92  double matchRatioPhi = 0;
93 
94  for (DTRecSegment4DCollection::const_iterator rechit = dtRecHits->begin(); rechit != dtRecHits->end(); ++rechit) {
95  LocalPoint pointLocal = rechit->localPosition();
96 
97  if (segments) {
98  // Loop over muon recHits
99  for (auto hit = dtHits.begin(); hit != dtHits.end(); ++hit) {
100  // Pick the one in the same DT Chamber as the muon
101  DetId idT = (*hit)->geographicalId();
102  if (!(rechit->geographicalId().rawId() == idT.rawId()))
103  continue;
104 
105  // and compare the local positions
106  LocalPoint segLocal = (*hit)->localPosition();
107  if ((fabs(pointLocal.x() - segLocal.x()) < ZCutParameter) &&
108  (fabs(pointLocal.y() - segLocal.y()) < ZCutParameter))
109  pointerTo4DSegments.push_back(&(*rechit));
110  }
111  continue;
112  }
113 
114  double nhitsPhi = 0;
115  double nhitsZ = 0;
116 
117  if (rechit->hasZed()) {
118  double countMuonDTHits = 0;
119  double countAgreeingHits = 0;
120 
121  const DTRecSegment2D* segmZ;
122  segmZ = dynamic_cast<const DTRecSegment2D*>(rechit->zSegment());
123  nhitsZ = segmZ->recHits().size();
124 
125  const vector<DTRecHit1D> hits1d = segmZ->specificRecHits();
126  DTChamberId chamberSegIdT((segmZ->geographicalId()).rawId());
127 
128  // Loop over muon recHits
129  for (auto hit = dtHits.begin(); hit != dtHits.end(); ++hit) {
130  if (!(*hit)->isValid())
131  continue;
132 
133  DetId idT = (*hit)->geographicalId();
134  DTChamberId dtDetIdHitT(idT.rawId());
135  DTSuperLayerId dtDetLayerIdHitT(idT.rawId());
136 
137  LocalPoint pointLocal = (*hit)->localPosition();
138 
139  if ((chamberSegIdT == dtDetIdHitT) && (dtDetLayerIdHitT.superlayer() == 2))
140  countMuonDTHits++;
141 
142  for (vector<DTRecHit1D>::const_iterator hiti = hits1d.begin(); hiti != hits1d.end(); hiti++) {
143  if (!hiti->isValid())
144  continue;
145 
146  // Pick the one in the same DT Layer as the 1D hit
147  if (!(hiti->geographicalId().rawId() == idT.rawId()))
148  continue;
149 
150  // and compare the local positions
151  LocalPoint segLocal = hiti->localPosition();
152  // cout << "Zed Segment Point = "<<pointLocal<<" Muon Point = "<<segLocal<<" Dist: "
153  // << (fabs(pointLocal.x()-segLocal.x()))+(fabs(pointLocal.y()-segLocal.y()))<< endl;
154  if ((fabs(pointLocal.x() - segLocal.x()) < ZCutParameter) &&
155  (fabs(pointLocal.y() - segLocal.y()) < ZCutParameter))
156  countAgreeingHits++;
157  } //End Segment Hit Iteration
158  } //End Muon Hit Iteration
159 
160  matchRatioZ = countMuonDTHits == 0 ? 0 : countAgreeingHits / countMuonDTHits;
161  if (nhitsZ)
162  if (countAgreeingHits / nhitsZ > matchRatioZ)
163  matchRatioZ = countAgreeingHits / nhitsZ;
164  } //End HasZed Check
165 
166  if (rechit->hasPhi()) {
167  double countMuonDTHits = 0;
168  double countAgreeingHits = 0;
169 
170  //PREPARE PARAMETERS FOR SEGMENT DETECTOR GEOMETRY
171  const DTRecSegment2D* segmPhi;
172  segmPhi = dynamic_cast<const DTRecSegment2D*>(rechit->phiSegment());
173  nhitsPhi = segmPhi->recHits().size();
174 
175  const vector<DTRecHit1D> hits1d = segmPhi->specificRecHits();
176  DTChamberId chamberSegIdT((segmPhi->geographicalId()).rawId());
177 
178  // Loop over muon recHits
179  for (auto hit = dtHits.begin(); hit != dtHits.end(); ++hit) {
180  if (!(*hit)->isValid())
181  continue;
182 
183  DetId idT = (*hit)->geographicalId();
184  DTChamberId dtDetIdHitT(idT.rawId());
185  DTSuperLayerId dtDetLayerIdHitT(idT.rawId());
186 
187  LocalPoint pointLocal =
188  (*hit)
189  ->localPosition(); //Localposition is in DTLayer http://cmslxr.fnal.gov/lxr/source/DataFormats/DTRecHit/interface/DTRecHit1D.h
190 
191  if ((chamberSegIdT == dtDetIdHitT) &&
192  ((dtDetLayerIdHitT.superlayer() == 1) || (dtDetLayerIdHitT.superlayer() == 3)))
193  countMuonDTHits++;
194 
195  for (vector<DTRecHit1D>::const_iterator hiti = hits1d.begin(); hiti != hits1d.end(); hiti++) {
196  if (!hiti->isValid())
197  continue;
198 
199  // Pick the one in the same DT Layer as the 1D hit
200  if (!(hiti->geographicalId().rawId() == idT.rawId()))
201  continue;
202 
203  // and compare the local positions
204  LocalPoint segLocal = hiti->localPosition();
205  // cout << " Phi Segment Point = "<<pointLocal<<" Muon Point = "<<segLocal<<" Dist: "
206  // << (fabs(pointLocal.x()-segLocal.x()))+(fabs(pointLocal.y()-segLocal.y()))<< endl;
207 
208  if ((fabs(pointLocal.x() - segLocal.x()) < PhiCutParameter) &&
209  (fabs(pointLocal.y() - segLocal.y()) < PhiCutParameter))
210  countAgreeingHits++;
211  } // End Segment Hit Iteration
212  } // End Muon Hit Iteration
213 
214  matchRatioPhi = countMuonDTHits != 0 ? countAgreeingHits / countMuonDTHits : 0;
215  if (nhitsPhi)
216  if (countAgreeingHits / nhitsPhi > matchRatioPhi)
217  matchRatioPhi = countAgreeingHits / nhitsPhi;
218  } // End HasPhi Check
219  // DTChamberId chamberSegId2((rechit->geographicalId()).rawId());
220  if (dtTightMatch && nhitsPhi && nhitsZ) {
221  if ((matchRatioPhi > 0.9) && (matchRatioZ > 0.9)) {
222  // cout<<"Making a tight match in Chamber "<<chamberSegId2<<endl;
223  pointerTo4DSegments.push_back(&(*rechit));
224  }
225  } else {
226  if ((matchRatioPhi > 0.9 && nhitsPhi) || (matchRatioZ > 0.9 && nhitsZ)) {
227  // cout<<"Making a loose match in Chamber "<<chamberSegId2<<endl;
228  pointerTo4DSegments.push_back(&(*rechit));
229  }
230  }
231 
232  } //End DT Segment Iteration
233 
234  return pointerTo4DSegments;
235 }
236 
237 vector<const CSCSegment*> MuonSegmentMatcher::matchCSC(const reco::Track& muon, const edm::Event& event) {
238  using namespace edm;
239 
240  edm::Handle<CSCSegmentCollection> allSegmentsCSC;
241  event.getByToken(allSegmentsCSCToken, allSegmentsCSC);
242 
243  vector<const CSCSegment*> pointerToCSCSegments;
244 
245  double matchRatioCSC = 0;
246  int numCSC = 0;
247  double CSCXCut = 0.001;
248  double CSCYCut = 0.001;
249  double countMuonCSCHits = 0;
250 
251  for (CSCSegmentCollection::const_iterator segmentCSC = allSegmentsCSC->begin(); segmentCSC != allSegmentsCSC->end();
252  segmentCSC++) {
253  double CSCcountAgreeingHits = 0;
254 
255  if (!segmentCSC->isValid())
256  continue;
257 
258  numCSC++;
259  const vector<CSCRecHit2D>& CSCRechits2D = segmentCSC->specificRecHits();
260  countMuonCSCHits = 0;
261  CSCDetId myChamber((*segmentCSC).geographicalId().rawId());
262 
263  bool segments = false;
264 
265  for (auto const& hitC : muon.recHits()) {
266  if (!hitC->isValid())
267  continue;
268  if (hitC->geographicalId().det() != DetId::Muon)
269  continue;
270  if (hitC->geographicalId().subdetId() != MuonSubdetId::CSC)
271  continue;
272  if (!hitC->isValid())
273  continue;
274  if (hitC->recHits().size() > 1)
275  segments = true;
276 
277  //DETECTOR CONSTRUCTION
278  DetId id = hitC->geographicalId();
279  CSCDetId cscDetIdHit(id.rawId());
280 
281  if (segments) {
282  if (!(myChamber.rawId() == cscDetIdHit.rawId()))
283  continue;
284 
285  // and compare the local positions
286  LocalPoint positionLocalCSC = hitC->localPosition();
287  LocalPoint segLocalCSC = segmentCSC->localPosition();
288  if ((fabs(positionLocalCSC.x() - segLocalCSC.x()) < CSCXCut) &&
289  (fabs(positionLocalCSC.y() - segLocalCSC.y()) < CSCYCut))
290  pointerToCSCSegments.push_back(&(*segmentCSC));
291  continue;
292  }
293 
294  if (!(cscDetIdHit.ring() == myChamber.ring()))
295  continue;
296  if (!(cscDetIdHit.station() == myChamber.station()))
297  continue;
298  if (!(cscDetIdHit.endcap() == myChamber.endcap()))
299  continue;
300  if (!(cscDetIdHit.chamber() == myChamber.chamber()))
301  continue;
302 
303  countMuonCSCHits++;
304 
305  LocalPoint positionLocalCSC = hitC->localPosition();
306 
307  for (vector<CSCRecHit2D>::const_iterator hiti = CSCRechits2D.begin(); hiti != CSCRechits2D.end(); hiti++) {
308  if (!hiti->isValid())
309  continue;
310  CSCDetId cscDetId((hiti->geographicalId()).rawId());
311 
312  if (hitC->geographicalId().rawId() != (hiti->geographicalId()).rawId())
313  continue;
314 
315  LocalPoint segLocalCSC = hiti->localPosition();
316  // cout<<"Layer Id (MuonHit) = "<<cscDetIdHit<<" Muon Local Position (det frame) "<<positionLocalCSC <<endl;
317  // cout<<"Layer Id (CSCHit) = "<<cscDetId<<" Hit Local Position (det frame) "<<segLocalCSC <<endl;
318  if ((fabs(positionLocalCSC.x() - segLocalCSC.x()) < CSCXCut) &&
319  (fabs(positionLocalCSC.y() - segLocalCSC.y()) < CSCYCut)) {
320  CSCcountAgreeingHits++;
321  // cout << " Matched." << endl;
322  }
323  } //End 2D rechit iteration
324  } //End muon hit iteration
325 
326  matchRatioCSC = countMuonCSCHits == 0 ? 0 : CSCcountAgreeingHits / countMuonCSCHits;
327 
328  if ((matchRatioCSC > 0.9) && ((countMuonCSCHits > 1) || !cscTightMatch))
329  pointerToCSCSegments.push_back(&(*segmentCSC));
330 
331  } //End CSC Segment Iteration
332 
333  return pointerToCSCSegments;
334 }
335 
336 vector<const RPCRecHit*> MuonSegmentMatcher::matchRPC(const reco::Track& muon, const edm::Event& event) {
337  using namespace edm;
338 
340  event.getByToken(rpcRecHitsToken, rpcRecHits);
341 
342  vector<const RPCRecHit*> pointerToRPCRecHits;
343  double RPCCut = 0.001;
344 
345  for (RPCRecHitCollection::const_iterator hitRPC = rpcRecHits->begin(); hitRPC != rpcRecHits->end(); hitRPC++) {
346  if (!hitRPC->isValid())
347  continue;
348 
349  RPCDetId myChamber((*hitRPC).geographicalId().rawId());
350  LocalPoint posLocalRPC = hitRPC->localPosition();
351  bool matched = false;
352 
353  for (auto const& hitC : muon.recHits()) {
354  if (!hitC->isValid())
355  continue;
356  if (hitC->geographicalId().det() != DetId::Muon)
357  continue;
358  if (hitC->geographicalId().subdetId() != MuonSubdetId::RPC)
359  continue;
360  if (!hitC->isValid())
361  continue;
362 
363  //DETECTOR CONSTRUCTION
364  DetId id = hitC->geographicalId();
365  RPCDetId rpcDetIdHit(id.rawId());
366 
367  if (rpcDetIdHit != myChamber)
368  continue;
369  LocalPoint posLocalMuon = hitC->localPosition();
370 
371  // cout<<"Layer Id (MuonHit) = "<<rpcDetIdHit<<" Muon Local Position (det frame) "<<posLocalMuon <<endl;
372  // cout<<"Layer Id (RPCHit) = "<<myChamber<<" Hit Local Position (det frame) "<<posLocalRPC <<endl;
373  if ((fabs(posLocalMuon.x() - posLocalRPC.x()) < RPCCut)) {
374  matched = true;
375  break;
376  }
377  }
378 
379  if (matched)
380  pointerToRPCRecHits.push_back(&(*hitRPC));
381  }
382 
383  return pointerToRPCRecHits;
384 }
385 
386 //define this as a plug-in
387 //DEFINE_FWK_MODULE(MuonSegmentMatcher);
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< RPCRecHitCollection > rpcRecHitsToken
MuonSegmentMatcher(const edm::ParameterSet &, edm::ConsumesCollector &iC)
constructor with Parameter Set and MuonServiceProxy
std::vector< const CSCSegment * > matchCSC(const reco::Track &muon, const edm::Event &event)
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:171
edm::InputTag RPCHitTags_
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
edm::EDGetTokenT< CSCSegmentCollection > allSegmentsCSCToken
int chamber() const
Definition: CSCDetId.h:62
std::vector< const TrackingRecHit * > recHits() const override
Access to component RecHits (if any)
edm::InputTag DTSegmentTags_
edm::InputTag CSCSegmentTags_
std::vector< const RPCRecHit * > matchRPC(const reco::Track &muon, const edm::Event &event)
virtual ~MuonSegmentMatcher()
destructor
Definition: DetId.h:17
int station() const
Definition: CSCDetId.h:79
DetId geographicalId() const
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
static constexpr int RPC
Definition: MuonSubdetId.h:13
int endcap() const
Definition: CSCDetId.h:85
HLT enums.
edm::EDGetTokenT< DTRecSegment4DCollection > dtRecHitsToken
std::vector< const DTRecSegment4D * > matchDT(const reco::Track &muon, const edm::Event &event)
perform the matching
static constexpr int DT
Definition: MuonSubdetId.h:11
int ring() const
Definition: CSCDetId.h:68
static constexpr int CSC
Definition: MuonSubdetId.h:12
Definition: event.py:1