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