CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
MuonSegmentMatcher Class Reference

#include <MuonSegmentMatcher.h>

Public Member Functions

std::vector< const CSCSegment * > matchCSC (const reco::Track &muon, const edm::Event &event)
 
std::vector< const DTRecSegment4D * > matchDT (const reco::Track &muon, const edm::Event &event)
 perform the matching More...
 
std::vector< const RPCRecHit * > matchRPC (const reco::Track &muon, const edm::Event &event)
 
 MuonSegmentMatcher (const edm::ParameterSet &, edm::ConsumesCollector &iC)
 constructor with Parameter Set and MuonServiceProxy More...
 
virtual ~MuonSegmentMatcher ()
 destructor More...
 

Private Attributes

edm::EDGetTokenT< CSCSegmentCollectionallSegmentsCSCToken
 
edm::InputTag CSCSegmentTags_
 
bool cscTightMatch
 
double dtRadius_
 
edm::EDGetTokenT< DTRecSegment4DCollectiondtRecHitsToken
 
edm::InputTag DTSegmentTags_
 
bool dtTightMatch
 
edm::InputTag RPCHitTags_
 
edm::EDGetTokenT< RPCRecHitCollectionrpcRecHitsToken
 
const edm::EventtheEvent
 
const MuonServiceProxytheService
 
edm::InputTag TKtrackTags_
 
edm::InputTag trackTags_
 

Detailed Description

Definition at line 24 of file MuonSegmentMatcher.h.

Constructor & Destructor Documentation

MuonSegmentMatcher::MuonSegmentMatcher ( const edm::ParameterSet matchParameters,
edm::ConsumesCollector iC 
)

constructor with Parameter Set and MuonServiceProxy

Definition at line 43 of file MuonSegmentMatcher.cc.

References allSegmentsCSCToken, edm::ConsumesCollector::consumes(), CSCSegmentTags_, dtRecHitsToken, DTSegmentTags_, edm::ParameterSet::existsAs(), edm::ParameterSet::getParameter(), RPCHitTags_, and rpcRecHitsToken.

44  :
45  DTSegmentTags_(matchParameters.getParameter<edm::InputTag>("DTsegments")),
46  CSCSegmentTags_(matchParameters.getParameter<edm::InputTag>("CSCsegments")),
47  dtRadius_(matchParameters.getParameter<double>("DTradius")),
48  dtTightMatch(matchParameters.getParameter<bool>("TightMatchDT")),
49  cscTightMatch(matchParameters.getParameter<bool>("TightMatchCSC"))
50 {
51  if (matchParameters.existsAs<edm::InputTag>("RPChits")) {
52  RPCHitTags_=matchParameters.getParameter<edm::InputTag>("RPChits");
53  } else RPCHitTags_=edm::InputTag("hltRpcRecHits");
54 
58 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:161
edm::EDGetTokenT< RPCRecHitCollection > rpcRecHitsToken
edm::InputTag RPCHitTags_
edm::EDGetTokenT< CSCSegmentCollection > allSegmentsCSCToken
edm::InputTag DTSegmentTags_
edm::InputTag CSCSegmentTags_
edm::EDGetTokenT< DTRecSegment4DCollection > dtRecHitsToken
MuonSegmentMatcher::~MuonSegmentMatcher ( )
virtual

destructor

Definition at line 60 of file MuonSegmentMatcher.cc.

61 {
62 }

Member Function Documentation

vector< const CSCSegment * > MuonSegmentMatcher::matchCSC ( const reco::Track muon,
const edm::Event event 
)

Definition at line 233 of file MuonSegmentMatcher.cc.

References allSegmentsCSCToken, CSCDetId::chamber(), MuonSubdetId::CSC, cscTightMatch, CSCDetId::endcap(), DetId::Muon, DetId::rawId(), reco::Track::recHits(), CSCDetId::ring(), CSCDetId::station(), PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by CSCHaloAlgo::Calculate(), and CSCTimingExtractor::fillTiming().

234 {
235 
236  using namespace edm;
237 
238  edm::Handle<CSCSegmentCollection> allSegmentsCSC;
239  event.getByToken(allSegmentsCSCToken, allSegmentsCSC);
240 
241  vector<const CSCSegment*> pointerToCSCSegments;
242 
243  double matchRatioCSC=0;
244  int numCSC = 0;
245  double CSCXCut = 0.001;
246  double CSCYCut = 0.001;
247  double countMuonCSCHits = 0;
248 
249  for(CSCSegmentCollection::const_iterator segmentCSC = allSegmentsCSC->begin(); segmentCSC != allSegmentsCSC->end(); segmentCSC++) {
250  double CSCcountAgreeingHits=0;
251 
252  if ( !segmentCSC->isValid()) continue;
253 
254  numCSC++;
255  const vector<CSCRecHit2D>& CSCRechits2D = segmentCSC->specificRecHits();
256  countMuonCSCHits = 0;
257  CSCDetId myChamber((*segmentCSC).geographicalId().rawId());
258 
259  bool segments = false;
260 
261  for(auto const& hitC : muon.recHits()) {
262  if (!hitC->isValid()) continue;
263  if ( hitC->geographicalId().det() != DetId::Muon ) continue;
264  if ( hitC->geographicalId().subdetId() != MuonSubdetId::CSC ) continue;
265  if (!hitC->isValid()) continue;
266  if ( hitC->recHits().size()>1) segments = true;
267 
268  //DETECTOR CONSTRUCTION
269  DetId id = hitC->geographicalId();
270  CSCDetId cscDetIdHit(id.rawId());
271 
272  if (segments) {
273  if(!(myChamber.rawId()==cscDetIdHit.rawId())) continue;
274 
275  // and compare the local positions
276  LocalPoint positionLocalCSC = hitC->localPosition();
277  LocalPoint segLocalCSC = segmentCSC->localPosition();
278  if ((fabs(positionLocalCSC.x()-segLocalCSC.x())<CSCXCut) &&
279  (fabs(positionLocalCSC.y()-segLocalCSC.y())<CSCYCut))
280  pointerToCSCSegments.push_back(&(*segmentCSC));
281  continue;
282  }
283 
284  if(!(cscDetIdHit.ring()==myChamber.ring())) continue;
285  if(!(cscDetIdHit.station()==myChamber.station())) continue;
286  if(!(cscDetIdHit.endcap()==myChamber.endcap())) continue;
287  if(!(cscDetIdHit.chamber()==myChamber.chamber())) continue;
288 
289  countMuonCSCHits++;
290 
291  LocalPoint positionLocalCSC = hitC->localPosition();
292 
293  for (vector<CSCRecHit2D>::const_iterator hiti=CSCRechits2D.begin(); hiti!=CSCRechits2D.end(); hiti++) {
294 
295  if ( !hiti->isValid()) continue;
296  CSCDetId cscDetId((hiti->geographicalId()).rawId());
297 
298  if (hitC->geographicalId().rawId()!=(hiti->geographicalId()).rawId()) continue;
299 
300  LocalPoint segLocalCSC = hiti->localPosition();
301  // cout<<"Layer Id (MuonHit) = "<<cscDetIdHit<<" Muon Local Position (det frame) "<<positionLocalCSC <<endl;
302  // cout<<"Layer Id (CSCHit) = "<<cscDetId<<" Hit Local Position (det frame) "<<segLocalCSC <<endl;
303  if((fabs(positionLocalCSC.x()-segLocalCSC.x())<CSCXCut) &&
304  (fabs(positionLocalCSC.y()-segLocalCSC.y())<CSCYCut)) {
305  CSCcountAgreeingHits++;
306  // cout << " Matched." << endl;
307  }
308  }//End 2D rechit iteration
309  }//End muon hit iteration
310 
311  matchRatioCSC = countMuonCSCHits == 0 ? 0 : CSCcountAgreeingHits/countMuonCSCHits;
312 
313  if ((matchRatioCSC>0.9) && ((countMuonCSCHits>1) || !cscTightMatch)) pointerToCSCSegments.push_back(&(*segmentCSC));
314 
315  } //End CSC Segment Iteration
316 
317  return pointerToCSCSegments;
318 
319 }
T y() const
Definition: PV3DBase.h:63
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
static const int CSC
Definition: MuonSubdetId.h:13
edm::EDGetTokenT< CSCSegmentCollection > allSegmentsCSCToken
auto recHits() const
Access to reconstructed hits on the track.
Definition: Track.h:106
Definition: DetId.h:18
HLT enums.
T x() const
Definition: PV3DBase.h:62
vector< const DTRecSegment4D * > MuonSegmentMatcher::matchDT ( const reco::Track muon,
const edm::Event event 
)

perform the matching

Definition at line 65 of file MuonSegmentMatcher.cc.

References MuonSubdetId::DT, dtRadius_, dtRecHitsToken, dtTightMatch, TrackingRecHit::geographicalId(), DetId::Muon, DetId::rawId(), reco::Track::recHits(), DTRecSegment2D::recHits(), DTRecSegment2D::specificRecHits(), PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by DTTimingExtractor::fillTiming().

66 {
67  using namespace edm;
68 
70  event.getByToken(dtRecHitsToken, dtRecHits);
71 
72  vector<const DTRecSegment4D*> pointerTo4DSegments;
73 
74  std::vector<TrackingRecHit const *> dtHits;
75 
76  bool segments = false;
77 
78  // Loop and select DT recHits
79  for(auto const& hit : muon.recHits()) {
80  if (!hit->isValid()) continue;
81  if ( hit->geographicalId().det() != DetId::Muon ) continue;
82  if ( hit->geographicalId().subdetId() != MuonSubdetId::DT ) continue;
83  if (!hit->recHits().empty())
84  if ((*hit->recHits().begin())->recHits().size()>1) segments = true;
85  dtHits.push_back(hit);
86  }
87 
88  // cout << "Muon DT hits found: " << dtHits.size() << " segments " << segments << endl;
89 
90  double PhiCutParameter=dtRadius_;
91  double ZCutParameter=dtRadius_;
92  double matchRatioZ=0;
93  double matchRatioPhi=0;
94 
95  for (DTRecSegment4DCollection::const_iterator rechit = dtRecHits->begin(); rechit!=dtRecHits->end();++rechit) {
96 
97  LocalPoint pointLocal = rechit->localPosition();
98 
99  if (segments) {
100  // Loop over muon recHits
101  for(auto hit = dtHits.begin(); hit != dtHits.end(); ++hit) {
102 
103  // Pick the one in the same DT Chamber as the muon
104  DetId idT = (*hit)->geographicalId();
105  if(!(rechit->geographicalId().rawId()==idT.rawId())) continue;
106 
107  // and compare the local positions
108  LocalPoint segLocal = (*hit)->localPosition();
109  if ((fabs(pointLocal.x()-segLocal.x())<ZCutParameter) &&
110  (fabs(pointLocal.y()-segLocal.y())<ZCutParameter))
111  pointerTo4DSegments.push_back(&(*rechit));
112  }
113  continue;
114  }
115 
116  double nhitsPhi = 0;
117  double nhitsZ = 0;
118 
119  if(rechit->hasZed()) {
120  double countMuonDTHits = 0;
121  double countAgreeingHits=0;
122 
123  const DTRecSegment2D* segmZ;
124  segmZ = dynamic_cast<const DTRecSegment2D*>(rechit->zSegment());
125  nhitsZ = segmZ->recHits().size();
126 
127  const vector<DTRecHit1D> hits1d = segmZ->specificRecHits();
128  DTChamberId chamberSegIdT((segmZ->geographicalId()).rawId());
129 
130  // Loop over muon recHits
131  for(auto hit = dtHits.begin(); hit != dtHits.end(); ++hit) {
132 
133  if ( !(*hit)->isValid()) continue;
134 
135  DetId idT = (*hit)->geographicalId();
136  DTChamberId dtDetIdHitT(idT.rawId());
137  DTSuperLayerId dtDetLayerIdHitT(idT.rawId());
138 
139  LocalPoint pointLocal = (*hit)->localPosition();
140 
141  if ((chamberSegIdT==dtDetIdHitT) && (dtDetLayerIdHitT.superlayer()==2)) countMuonDTHits++;
142 
143  for (vector<DTRecHit1D>::const_iterator hiti=hits1d.begin(); hiti!=hits1d.end(); hiti++) {
144 
145  if ( !hiti->isValid()) continue;
146 
147  // Pick the one in the same DT Layer as the 1D hit
148  if(!(hiti->geographicalId().rawId()==idT.rawId())) 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) 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 
180  if ( !(*hit)->isValid()) continue;
181 
182  DetId idT = (*hit)->geographicalId();
183  DTChamberId dtDetIdHitT(idT.rawId());
184  DTSuperLayerId dtDetLayerIdHitT(idT.rawId());
185 
186  LocalPoint pointLocal = (*hit)->localPosition(); //Localposition is in DTLayer http://cmslxr.fnal.gov/lxr/source/DataFormats/DTRecHit/interface/DTRecHit1D.h
187 
188  if ((chamberSegIdT==dtDetIdHitT)&&((dtDetLayerIdHitT.superlayer()==1)||(dtDetLayerIdHitT.superlayer()==3)))
189  countMuonDTHits++;
190 
191  for (vector<DTRecHit1D>::const_iterator hiti=hits1d.begin(); hiti!=hits1d.end(); hiti++) {
192 
193  if ( !hiti->isValid()) continue;
194 
195  // Pick the one in the same DT Layer as the 1D hit
196  if(!(hiti->geographicalId().rawId()==idT.rawId())) continue;
197 
198  // and compare the local positions
199  LocalPoint segLocal = hiti->localPosition();
200 // cout << " Phi Segment Point = "<<pointLocal<<" Muon Point = "<<segLocal<<" Dist: "
201 // << (fabs(pointLocal.x()-segLocal.x()))+(fabs(pointLocal.y()-segLocal.y()))<< endl;
202 
203  if ((fabs(pointLocal.x()-segLocal.x())<PhiCutParameter) &&
204  (fabs(pointLocal.y()-segLocal.y())<PhiCutParameter))
205  countAgreeingHits++;
206  } // End Segment Hit Iteration
207  } // End Muon Hit Iteration
208 
209  matchRatioPhi = countMuonDTHits != 0 ? countAgreeingHits/countMuonDTHits : 0;
210  if (nhitsPhi)
211  if (countAgreeingHits/nhitsPhi>matchRatioPhi) matchRatioPhi=countAgreeingHits/nhitsPhi;
212  } // End HasPhi Check
213 // DTChamberId chamberSegId2((rechit->geographicalId()).rawId());
214  if (dtTightMatch && nhitsPhi && nhitsZ) {
215  if((matchRatioPhi>0.9)&&(matchRatioZ>0.9)) {
216 // cout<<"Making a tight match in Chamber "<<chamberSegId2<<endl;
217  pointerTo4DSegments.push_back(&(*rechit));
218  }
219  } else {
220  if((matchRatioPhi>0.9 && nhitsPhi)||(matchRatioZ>0.9 && nhitsZ)) {
221 // cout<<"Making a loose match in Chamber "<<chamberSegId2<<endl;
222  pointerTo4DSegments.push_back(&(*rechit));
223  }
224  }
225 
226  } //End DT Segment Iteration
227 
228  return pointerTo4DSegments;
229 }
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
T y() const
Definition: PV3DBase.h:63
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
auto recHits() const
Access to reconstructed hits on the track.
Definition: Track.h:106
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
Definition: DetId.h:18
HLT enums.
edm::EDGetTokenT< DTRecSegment4DCollection > dtRecHitsToken
std::vector< const TrackingRecHit * > recHits() const override
Access to component RecHits (if any)
static const int DT
Definition: MuonSubdetId.h:12
DetId geographicalId() const
T x() const
Definition: PV3DBase.h:62
vector< const RPCRecHit * > MuonSegmentMatcher::matchRPC ( const reco::Track muon,
const edm::Event event 
)

Definition at line 322 of file MuonSegmentMatcher.cc.

References electrons_cff::matched, DetId::Muon, reco::Track::recHits(), MuonSubdetId::RPC, JetIDParams_cfi::rpcRecHits, rpcRecHitsToken, and PV3DBase< T, PVType, FrameType >::x().

323 {
324 
325  using namespace edm;
326 
328  event.getByToken(rpcRecHitsToken, rpcRecHits);
329 
330  vector<const RPCRecHit*> pointerToRPCRecHits;
331  double RPCCut = 0.001;
332 
333  for(RPCRecHitCollection::const_iterator hitRPC = rpcRecHits->begin(); hitRPC != rpcRecHits->end(); hitRPC++) {
334 
335  if ( !hitRPC->isValid()) continue;
336 
337  RPCDetId myChamber((*hitRPC).geographicalId().rawId());
338  LocalPoint posLocalRPC = hitRPC->localPosition();
339  bool matched=false;
340 
341  for(auto const& hitC : muon.recHits()) {
342  if (!hitC->isValid()) continue;
343  if ( hitC->geographicalId().det() != DetId::Muon ) continue;
344  if ( hitC->geographicalId().subdetId() != MuonSubdetId::RPC ) continue;
345  if (!hitC->isValid()) continue;
346 
347  //DETECTOR CONSTRUCTION
348  DetId id = hitC->geographicalId();
349  RPCDetId rpcDetIdHit(id.rawId());
350 
351  if (rpcDetIdHit!=myChamber) continue;
352  LocalPoint posLocalMuon = hitC->localPosition();
353 
354 // cout<<"Layer Id (MuonHit) = "<<rpcDetIdHit<<" Muon Local Position (det frame) "<<posLocalMuon <<endl;
355 // cout<<"Layer Id (RPCHit) = "<<myChamber<<" Hit Local Position (det frame) "<<posLocalRPC <<endl;
356  if((fabs(posLocalMuon.x()-posLocalRPC.x())<RPCCut)) {
357  matched=true;
358  break;
359  }
360 
361  }
362 
363  if (matched) pointerToRPCRecHits.push_back(&(*hitRPC));
364  }
365 
366  return pointerToRPCRecHits;
367 }
edm::EDGetTokenT< RPCRecHitCollection > rpcRecHitsToken
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
auto recHits() const
Access to reconstructed hits on the track.
Definition: Track.h:106
Definition: DetId.h:18
HLT enums.
static const int RPC
Definition: MuonSubdetId.h:14
T x() const
Definition: PV3DBase.h:62

Member Data Documentation

edm::EDGetTokenT<CSCSegmentCollection> MuonSegmentMatcher::allSegmentsCSCToken
private

Definition at line 54 of file MuonSegmentMatcher.h.

Referenced by matchCSC(), and MuonSegmentMatcher().

edm::InputTag MuonSegmentMatcher::CSCSegmentTags_
private

Definition at line 50 of file MuonSegmentMatcher.h.

Referenced by MuonSegmentMatcher().

bool MuonSegmentMatcher::cscTightMatch
private

Definition at line 60 of file MuonSegmentMatcher.h.

Referenced by matchCSC().

double MuonSegmentMatcher::dtRadius_
private

Definition at line 57 of file MuonSegmentMatcher.h.

Referenced by matchDT().

edm::EDGetTokenT<DTRecSegment4DCollection> MuonSegmentMatcher::dtRecHitsToken
private

Definition at line 53 of file MuonSegmentMatcher.h.

Referenced by matchDT(), and MuonSegmentMatcher().

edm::InputTag MuonSegmentMatcher::DTSegmentTags_
private

Definition at line 49 of file MuonSegmentMatcher.h.

Referenced by MuonSegmentMatcher().

bool MuonSegmentMatcher::dtTightMatch
private

Definition at line 59 of file MuonSegmentMatcher.h.

Referenced by matchDT().

edm::InputTag MuonSegmentMatcher::RPCHitTags_
private

Definition at line 51 of file MuonSegmentMatcher.h.

Referenced by MuonSegmentMatcher().

edm::EDGetTokenT<RPCRecHitCollection> MuonSegmentMatcher::rpcRecHitsToken
private

Definition at line 55 of file MuonSegmentMatcher.h.

Referenced by matchRPC(), and MuonSegmentMatcher().

const edm::Event* MuonSegmentMatcher::theEvent
private

Definition at line 45 of file MuonSegmentMatcher.h.

const MuonServiceProxy* MuonSegmentMatcher::theService
private

Definition at line 44 of file MuonSegmentMatcher.h.

edm::InputTag MuonSegmentMatcher::TKtrackTags_
private

Definition at line 47 of file MuonSegmentMatcher.h.

edm::InputTag MuonSegmentMatcher::trackTags_
private

Definition at line 48 of file MuonSegmentMatcher.h.