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 *, edm::ConsumesCollector &iC)
 constructor with Parameter Set and MuonServiceProxy More...
 
virtual ~MuonSegmentMatcher ()
 destructor More...
 

Private Attributes

edm::EDGetTokenT
< CSCSegmentCollection
allSegmentsCSCToken
 
edm::InputTag CSCSegmentTags_
 
bool cscTightMatch
 
double dtRadius_
 
edm::EDGetTokenT
< DTRecSegment4DCollection
dtRecHitsToken
 
edm::InputTag DTSegmentTags_
 
bool dtTightMatch
 
const edm::EventtheEvent
 
const MuonServiceProxytheService
 
edm::InputTag TKtrackTags_
 
edm::InputTag trackTags_
 

Detailed Description

Definition at line 23 of file MuonSegmentMatcher.h.

Constructor & Destructor Documentation

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

constructor with Parameter Set and MuonServiceProxy

Definition at line 45 of file MuonSegmentMatcher.cc.

References allSegmentsCSCToken, edm::ConsumesCollector::consumes(), CSCSegmentTags_, dtRecHitsToken, and DTSegmentTags_.

46  :
47  theService(service),
48  DTSegmentTags_(matchParameters.getParameter<edm::InputTag>("DTsegments")),
49  CSCSegmentTags_(matchParameters.getParameter<edm::InputTag>("CSCsegments")),
50  dtRadius_(matchParameters.getParameter<double>("DTradius")),
51  dtTightMatch(matchParameters.getParameter<bool>("TightMatchDT")),
52  cscTightMatch(matchParameters.getParameter<bool>("TightMatchCSC"))
53 {
56 
57 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
const MuonServiceProxy * theService
edm::EDGetTokenT< CSCSegmentCollection > allSegmentsCSCToken
edm::InputTag DTSegmentTags_
edm::InputTag CSCSegmentTags_
edm::EDGetTokenT< DTRecSegment4DCollection > dtRecHitsToken
MuonSegmentMatcher::~MuonSegmentMatcher ( )
virtual

destructor

Definition at line 59 of file MuonSegmentMatcher.cc.

60 {
61 }

Member Function Documentation

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

Definition at line 232 of file MuonSegmentMatcher.cc.

References allSegmentsCSCToken, CSCDetId::chamber(), MuonSubdetId::CSC, 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().

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

References MuonSubdetId::DT, dtRadius_, dtRecHitsToken, 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().

65 {
66  using namespace edm;
67 
69  event.getByToken(dtRecHitsToken, dtRecHits);
70 
71  vector<const DTRecSegment4D*> pointerTo4DSegments;
72 
74 
75  bool segments = false;
76 
77  // Loop and select DT recHits
78  for(trackingRecHit_iterator hit = muon.recHitsBegin(); hit != muon.recHitsEnd(); ++hit) {
79  if (!(*hit)->isValid()) continue;
80  if ( (*hit)->geographicalId().det() != DetId::Muon ) continue;
81  if ( (*hit)->geographicalId().subdetId() != MuonSubdetId::DT ) continue;
82  if ((*hit)->recHits().size())
83  if ((*(*hit)->recHits().begin())->recHits().size()>1) 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 
96  LocalPoint pointLocal = rechit->localPosition();
97 
98  if (segments) {
99  // Loop over muon recHits
100  for(trackingRecHit_iterator hit = dtHits.begin(); hit != dtHits.end(); ++hit) {
101 
102  // Pick the one in the same DT Chamber as the muon
103  DetId idT = (*hit)->geographicalId();
104  if(!(rechit->geographicalId().rawId()==idT.rawId())) continue;
105 
106  // and compare the local positions
107  LocalPoint segLocal = (*hit)->localPosition();
108  if ((fabs(pointLocal.x()-segLocal.x())<ZCutParameter) &&
109  (fabs(pointLocal.y()-segLocal.y())<ZCutParameter))
110  pointerTo4DSegments.push_back(&(*rechit));
111  }
112  continue;
113  }
114 
115  double nhitsPhi = 0;
116  double nhitsZ = 0;
117 
118  if(rechit->hasZed()) {
119  double countMuonDTHits = 0;
120  double countAgreeingHits=0;
121 
122  const DTRecSegment2D* segmZ;
123  segmZ = dynamic_cast<const DTRecSegment2D*>(rechit->zSegment());
124  nhitsZ = segmZ->recHits().size();
125 
126  const vector<DTRecHit1D> hits1d = segmZ->specificRecHits();
127  DTChamberId chamberSegIdT((segmZ->geographicalId()).rawId());
128 
129  // Loop over muon recHits
130  for(trackingRecHit_iterator hit = dtHits.begin(); hit != dtHits.end(); ++hit) {
131 
132  if ( !(*hit)->isValid()) continue;
133 
134  DetId idT = (*hit)->geographicalId();
135  DTChamberId dtDetIdHitT(idT.rawId());
136  DTSuperLayerId dtDetLayerIdHitT(idT.rawId());
137 
138  LocalPoint pointLocal = (*hit)->localPosition();
139 
140  if ((chamberSegIdT==dtDetIdHitT) && (dtDetLayerIdHitT.superlayer()==2)) countMuonDTHits++;
141 
142  for (vector<DTRecHit1D>::const_iterator hiti=hits1d.begin(); hiti!=hits1d.end(); hiti++) {
143 
144  if ( !hiti->isValid()) continue;
145 
146  // Pick the one in the same DT Layer as the 1D hit
147  if(!(hiti->geographicalId().rawId()==idT.rawId())) 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) matchRatioZ=countAgreeingHits/nhitsZ;
162  } //End HasZed Check
163 
164  if(rechit->hasPhi()) {
165  double countMuonDTHits = 0;
166  double countAgreeingHits=0;
167 
168  //PREPARE PARAMETERS FOR SEGMENT DETECTOR GEOMETRY
169  const DTRecSegment2D* segmPhi;
170  segmPhi = dynamic_cast<const DTRecSegment2D*>(rechit->phiSegment());
171  nhitsPhi = segmPhi->recHits().size();
172 
173  const vector<DTRecHit1D> hits1d = segmPhi->specificRecHits();
174  DTChamberId chamberSegIdT((segmPhi->geographicalId()).rawId());
175 
176  // Loop over muon recHits
177  for(trackingRecHit_iterator hit = dtHits.begin(); hit != dtHits.end(); ++hit) {
178 
179  if ( !(*hit)->isValid()) continue;
180 
181  DetId idT = (*hit)->geographicalId();
182  DTChamberId dtDetIdHitT(idT.rawId());
183  DTSuperLayerId dtDetLayerIdHitT(idT.rawId());
184 
185  LocalPoint pointLocal = (*hit)->localPosition(); //Localposition is in DTLayer http://cmslxr.fnal.gov/lxr/source/DataFormats/DTRecHit/interface/DTRecHit1D.h
186 
187  if ((chamberSegIdT==dtDetIdHitT)&&((dtDetLayerIdHitT.superlayer()==1)||(dtDetLayerIdHitT.superlayer()==3)))
188  countMuonDTHits++;
189 
190  for (vector<DTRecHit1D>::const_iterator hiti=hits1d.begin(); hiti!=hits1d.end(); hiti++) {
191 
192  if ( !hiti->isValid()) continue;
193 
194  // Pick the one in the same DT Layer as the 1D hit
195  if(!(hiti->geographicalId().rawId()==idT.rawId())) continue;
196 
197  // and compare the local positions
198  LocalPoint segLocal = hiti->localPosition();
199 // cout << " Phi Segment Point = "<<pointLocal<<" Muon Point = "<<segLocal<<" Dist: "
200 // << (fabs(pointLocal.x()-segLocal.x()))+(fabs(pointLocal.y()-segLocal.y()))<< endl;
201 
202  if ((fabs(pointLocal.x()-segLocal.x())<PhiCutParameter) &&
203  (fabs(pointLocal.y()-segLocal.y())<PhiCutParameter))
204  countAgreeingHits++;
205  } // End Segment Hit Iteration
206  } // End Muon Hit Iteration
207 
208  matchRatioPhi = countMuonDTHits != 0 ? countAgreeingHits/countMuonDTHits : 0;
209  if (nhitsPhi)
210  if (countAgreeingHits/nhitsPhi>matchRatioPhi) matchRatioPhi=countAgreeingHits/nhitsPhi;
211  } // End HasPhi Check
212 // DTChamberId chamberSegId2((rechit->geographicalId()).rawId());
213  if (dtTightMatch && nhitsPhi && nhitsZ) {
214  if((matchRatioPhi>0.9)&&(matchRatioZ>0.9)) {
215 // cout<<"Making a tight match in Chamber "<<chamberSegId2<<endl;
216  pointerTo4DSegments.push_back(&(*rechit));
217  }
218  } else {
219  if((matchRatioPhi>0.9 && nhitsPhi)||(matchRatioZ>0.9 && nhitsZ)) {
220 // cout<<"Making a loose match in Chamber "<<chamberSegId2<<endl;
221  pointerTo4DSegments.push_back(&(*rechit));
222  }
223  }
224 
225  } //End DT Segment Iteration
226 
227  return pointerTo4DSegments;
228 }
T y() const
Definition: PV3DBase.h:63
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
virtual std::vector< const TrackingRecHit * > recHits() const
Access to component RecHits (if any)
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:18
edm::EDGetTokenT< DTRecSegment4DCollection > dtRecHitsToken
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:12
DetId geographicalId() const
T x() const
Definition: PV3DBase.h:62
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:65

Member Data Documentation

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

Definition at line 50 of file MuonSegmentMatcher.h.

Referenced by matchCSC(), and MuonSegmentMatcher().

edm::InputTag MuonSegmentMatcher::CSCSegmentTags_
private

Definition at line 47 of file MuonSegmentMatcher.h.

Referenced by MuonSegmentMatcher().

bool MuonSegmentMatcher::cscTightMatch
private

Definition at line 56 of file MuonSegmentMatcher.h.

Referenced by matchCSC().

double MuonSegmentMatcher::dtRadius_
private

Definition at line 53 of file MuonSegmentMatcher.h.

Referenced by matchDT().

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

Definition at line 49 of file MuonSegmentMatcher.h.

Referenced by matchDT(), and MuonSegmentMatcher().

edm::InputTag MuonSegmentMatcher::DTSegmentTags_
private

Definition at line 46 of file MuonSegmentMatcher.h.

Referenced by MuonSegmentMatcher().

bool MuonSegmentMatcher::dtTightMatch
private

Definition at line 55 of file MuonSegmentMatcher.h.

Referenced by matchDT().

const edm::Event* MuonSegmentMatcher::theEvent
private

Definition at line 42 of file MuonSegmentMatcher.h.

const MuonServiceProxy* MuonSegmentMatcher::theService
private

Definition at line 41 of file MuonSegmentMatcher.h.

edm::InputTag MuonSegmentMatcher::TKtrackTags_
private

Definition at line 44 of file MuonSegmentMatcher.h.

edm::InputTag MuonSegmentMatcher::trackTags_
private

Definition at line 45 of file MuonSegmentMatcher.h.