CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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...
 
 MuonSegmentMatcher (const edm::ParameterSet &, MuonServiceProxy *)
 constructor with Parameter Set and MuonServiceProxy More...
 
virtual ~MuonSegmentMatcher ()
 destructor More...
 

Private Attributes

edm::InputTag CSCSegmentTags_
 
bool cscTightMatch
 
double dtRadius_
 
edm::InputTag DTSegmentTags_
 
bool dtTightMatch
 
const edm::EventtheEvent
 
const MuonServiceProxytheService
 
edm::InputTag TKtrackTags_
 
edm::InputTag trackTags_
 

Detailed Description

Definition at line 22 of file MuonSegmentMatcher.h.

Constructor & Destructor Documentation

MuonSegmentMatcher::MuonSegmentMatcher ( const edm::ParameterSet matchParameters,
MuonServiceProxy service 
)

constructor with Parameter Set and MuonServiceProxy

Definition at line 46 of file MuonSegmentMatcher.cc.

47  :
48  theService(service),
49  DTSegmentTags_(matchParameters.getParameter<edm::InputTag>("DTsegments")),
50  CSCSegmentTags_(matchParameters.getParameter<edm::InputTag>("CSCsegments")),
51  dtRadius_(matchParameters.getParameter<double>("DTradius")),
52  dtTightMatch(matchParameters.getParameter<bool>("TightMatchDT")),
53  cscTightMatch(matchParameters.getParameter<bool>("TightMatchCSC"))
54 {
55 }
T getParameter(std::string const &) const
const MuonServiceProxy * theService
edm::InputTag DTSegmentTags_
edm::InputTag CSCSegmentTags_
MuonSegmentMatcher::~MuonSegmentMatcher ( )
virtual

destructor

Definition at line 57 of file MuonSegmentMatcher.cc.

58 {
59 }

Member Function Documentation

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

Definition at line 230 of file MuonSegmentMatcher.cc.

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

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

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

perform the matching

Definition at line 62 of file MuonSegmentMatcher.cc.

References MuonSubdetId::DT, dtRadius_, DTSegmentTags_, dtTightMatch, TrackingRecHit::geographicalId(), DetId::Muon, edm::RefVector< C, T, F >::push_back(), DetId::rawId(), DTRecSegment2D::recHits(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), DTRecSegment2D::specificRecHits(), PV3DBase< T, PVType, FrameType >::x(), and PV3DBase< T, PVType, FrameType >::y().

Referenced by DTTimingExtractor::fillTiming().

63 {
64  using namespace edm;
65 
67  event.getByLabel(DTSegmentTags_, dtRecHits);
68 
69  vector<const DTRecSegment4D*> pointerTo4DSegments;
70 
72 
73  bool segments = false;
74 
75  // Loop and select DT recHits
76  for(trackingRecHit_iterator hit = muon.recHitsBegin(); hit != muon.recHitsEnd(); ++hit) {
77  if ( !(*hit)->isValid()) continue;
78  if ( (*hit)->geographicalId().det() != DetId::Muon ) continue;
79  if ( (*hit)->geographicalId().subdetId() != MuonSubdetId::DT ) continue;
80  if (!(*hit)->isValid()) continue;
81  if ((*hit)->recHits().size()>1) segments = true;
82  dtHits.push_back(*hit);
83  }
84 
85  double PhiCutParameter=dtRadius_;
86  double ZCutParameter=dtRadius_;
87  double matchRatioZ=0;
88  double matchRatioPhi=0;
89 
90  for (DTRecSegment4DCollection::const_iterator rechit = dtRecHits->begin(); rechit!=dtRecHits->end();++rechit) {
91 
92  if ( !rechit->isValid()) continue;
93  LocalPoint pointLocal = rechit->localPosition();
94 
95  if (segments) {
96  // Loop over muon recHits
97  for(trackingRecHit_iterator hit = dtHits.begin(); hit != dtHits.end(); ++hit) {
98  if ( !(*hit)->isValid()) continue;
99 
100  // Pick the one in the same DT Chamber as the muon
101  DetId idT = (*hit)->geographicalId();
102  if(!(rechit->geographicalId().rawId()==idT.rawId())) 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(trackingRecHit_iterator hit = dtHits.begin(); hit != dtHits.end(); ++hit) {
129 
130  if ( !(*hit)->isValid()) 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)) countMuonDTHits++;
139 
140  for (vector<DTRecHit1D>::const_iterator hiti=hits1d.begin(); hiti!=hits1d.end(); hiti++) {
141 
142  if ( !hiti->isValid()) continue;
143 
144  // Pick the one in the same DT Layer as the 1D hit
145  if(!(hiti->geographicalId().rawId()==idT.rawId())) continue;
146 
147  // and compare the local positions
148  LocalPoint segLocal = hiti->localPosition();
149 // cout << "Zed Segment Point = "<<pointLocal<<" Muon Point = "<<segLocal<<" Dist: "
150 // << (fabs(pointLocal.x()-segLocal.x()))+(fabs(pointLocal.y()-segLocal.y()))<< endl;
151  if ((fabs(pointLocal.x()-segLocal.x())<ZCutParameter) &&
152  (fabs(pointLocal.y()-segLocal.y())<ZCutParameter))
153  countAgreeingHits++;
154  } //End Segment Hit Iteration
155  } //End Muon Hit Iteration
156 
157  matchRatioZ = countMuonDTHits == 0 ? 0 : countAgreeingHits/countMuonDTHits;
158  if (nhitsZ)
159  if (countAgreeingHits/nhitsZ>matchRatioZ) matchRatioZ=countAgreeingHits/nhitsZ;
160  } //End HasZed Check
161 
162  if(rechit->hasPhi()) {
163  double countMuonDTHits = 0;
164  double countAgreeingHits=0;
165 
166  //PREPARE PARAMETERS FOR SEGMENT DETECTOR GEOMETRY
167  const DTRecSegment2D* segmPhi;
168  segmPhi = dynamic_cast<const DTRecSegment2D*>(rechit->phiSegment());
169  nhitsPhi = segmPhi->recHits().size();
170 
171  const vector<DTRecHit1D> hits1d = segmPhi->specificRecHits();
172  DTChamberId chamberSegIdT((segmPhi->geographicalId()).rawId());
173 
174  // Loop over muon recHits
175  for(trackingRecHit_iterator hit = dtHits.begin(); hit != dtHits.end(); ++hit) {
176 
177  if ( !(*hit)->isValid()) continue;
178 
179  DetId idT = (*hit)->geographicalId();
180  DTChamberId dtDetIdHitT(idT.rawId());
181  DTSuperLayerId dtDetLayerIdHitT(idT.rawId());
182 
183  LocalPoint pointLocal = (*hit)->localPosition(); //Localposition is in DTLayer http://cmslxr.fnal.gov/lxr/source/DataFormats/DTRecHit/interface/DTRecHit1D.h
184 
185  if ((chamberSegIdT==dtDetIdHitT)&&((dtDetLayerIdHitT.superlayer()==1)||(dtDetLayerIdHitT.superlayer()==3)))
186  countMuonDTHits++;
187 
188  for (vector<DTRecHit1D>::const_iterator hiti=hits1d.begin(); hiti!=hits1d.end(); hiti++) {
189 
190  if ( !hiti->isValid()) continue;
191 
192  // Pick the one in the same DT Layer as the 1D hit
193  if(!(hiti->geographicalId().rawId()==idT.rawId())) continue;
194 
195  // and compare the local positions
196  LocalPoint segLocal = hiti->localPosition();
197 // cout << " Phi Segment Point = "<<pointLocal<<" Muon Point = "<<segLocal<<" Dist: "
198 // << (fabs(pointLocal.x()-segLocal.x()))+(fabs(pointLocal.y()-segLocal.y()))<< endl;
199 
200  if ((fabs(pointLocal.x()-segLocal.x())<PhiCutParameter) &&
201  (fabs(pointLocal.y()-segLocal.y())<PhiCutParameter))
202  countAgreeingHits++;
203  } // End Segment Hit Iteration
204  } // End Muon Hit Iteration
205 
206  matchRatioPhi = countMuonDTHits != 0 ? countAgreeingHits/countMuonDTHits : 0;
207  if (nhitsPhi)
208  if (countAgreeingHits/nhitsPhi>matchRatioPhi) matchRatioPhi=countAgreeingHits/nhitsPhi;
209  } // End HasPhi Check
210 // DTChamberId chamberSegId2((rechit->geographicalId()).rawId());
211  if (dtTightMatch && nhitsPhi && nhitsZ) {
212  if((matchRatioPhi>0.9)&&(matchRatioZ>0.9)) {
213 // cout<<"Making a tight match in Chamber "<<chamberSegId2<<endl;
214  pointerTo4DSegments.push_back(&(*rechit));
215  }
216  } else {
217  if((matchRatioPhi>0.9 && nhitsPhi)||(matchRatioZ>0.9 && nhitsZ)) {
218 // cout<<"Making a loose match in Chamber "<<chamberSegId2<<endl;
219  pointerTo4DSegments.push_back(&(*rechit));
220  }
221  }
222 
223  } //End DT Segment Iteration
224 
225  return pointerTo4DSegments;
226 }
T y() const
Definition: PV3DBase.h:62
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:45
virtual std::vector< const TrackingRecHit * > recHits() const
Access to component RecHits (if any)
edm::InputTag DTSegmentTags_
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
Definition: DetId.h:20
void push_back(value_type const &ref)
Add a Ref&lt;C, T&gt; to the RefVector.
Definition: RefVector.h:64
static const int DT
Definition: MuonSubdetId.h:14
DetId geographicalId() const
T x() const
Definition: PV3DBase.h:61
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:65

Member Data Documentation

edm::InputTag MuonSegmentMatcher::CSCSegmentTags_
private

Definition at line 46 of file MuonSegmentMatcher.h.

Referenced by matchCSC().

bool MuonSegmentMatcher::cscTightMatch
private

Definition at line 50 of file MuonSegmentMatcher.h.

Referenced by matchCSC().

double MuonSegmentMatcher::dtRadius_
private

Definition at line 47 of file MuonSegmentMatcher.h.

Referenced by matchDT().

edm::InputTag MuonSegmentMatcher::DTSegmentTags_
private

Definition at line 45 of file MuonSegmentMatcher.h.

Referenced by matchDT().

bool MuonSegmentMatcher::dtTightMatch
private

Definition at line 49 of file MuonSegmentMatcher.h.

Referenced by matchDT().

const edm::Event* MuonSegmentMatcher::theEvent
private

Definition at line 41 of file MuonSegmentMatcher.h.

const MuonServiceProxy* MuonSegmentMatcher::theService
private

Definition at line 40 of file MuonSegmentMatcher.h.

edm::InputTag MuonSegmentMatcher::TKtrackTags_
private

Definition at line 43 of file MuonSegmentMatcher.h.

edm::InputTag MuonSegmentMatcher::trackTags_
private

Definition at line 44 of file MuonSegmentMatcher.h.