Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013 #include "RecoMuon/TrackingTools/interface/MuonSegmentMatcher.h"
00014
00015
00016
00017 #include "FWCore/Framework/interface/EDAnalyzer.h"
00018 #include "FWCore/Framework/interface/Event.h"
00019
00020 #include "DataFormats/MuonReco/interface/Muon.h"
00021 #include "RecoMuon/TrackingTools/interface/MuonServiceProxy.h"
00022
00023 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00024 #include "DataFormats/TrackReco/interface/Track.h"
00025
00026 #include "DataFormats/MuonDetId/interface/MuonSubdetId.h"
00027 #include "DataFormats/MuonDetId/interface/CSCDetId.h"
00028
00029 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
00030 #include "DataFormats/DTRecHit/interface/DTRecHitCollection.h"
00031 #include "DataFormats/CSCRecHit/interface/CSCSegmentCollection.h"
00032 #include "DataFormats/CSCRecHit/interface/CSCRecHit2DCollection.h"
00033 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00034
00035
00036 #include <vector>
00037 #include <iostream>
00038
00039 class MuonServiceProxy;
00040 class MuonSegmentMatcher;
00041
00042 using namespace std;
00043
00044
00045
00046 MuonSegmentMatcher::MuonSegmentMatcher(const edm::ParameterSet& matchParameters, MuonServiceProxy* service)
00047 :
00048 theService(service),
00049 DTSegmentTags_(matchParameters.getParameter<edm::InputTag>("DTsegments")),
00050 CSCSegmentTags_(matchParameters.getParameter<edm::InputTag>("CSCsegments")),
00051 dtTightMatch(matchParameters.getParameter<bool>("TightMatchDT")),
00052 cscTightMatch(matchParameters.getParameter<bool>("TightMatchCSC"))
00053 {
00054 }
00055
00056 MuonSegmentMatcher::~MuonSegmentMatcher()
00057 {
00058 }
00059
00060
00061 vector<const DTRecSegment4D*> MuonSegmentMatcher::matchDT(const reco::Track &muon, const edm::Event& event)
00062 {
00063 using namespace edm;
00064
00065 edm::Handle<DTRecSegment4DCollection> dtRecHits;
00066 event.getByLabel(DTSegmentTags_, dtRecHits);
00067
00068 vector<const DTRecSegment4D*> pointerTo4DSegments;
00069
00070 TrackingRecHitRefVector dtHits;
00071
00072
00073 for(trackingRecHit_iterator hit = muon.recHitsBegin(); hit != muon.recHitsEnd(); ++hit) {
00074 if ( !(*hit)->isValid()) continue;
00075 if ( (*hit)->geographicalId().det() != DetId::Muon ) continue;
00076 if ( (*hit)->geographicalId().subdetId() != MuonSubdetId::DT ) continue;
00077 if (!(*hit)->isValid()) continue;
00078 dtHits.push_back(*hit);
00079 }
00080
00081 double PhiCutParameter=0.001;
00082 double ZCutParameter=0.001;
00083 double matchRatioZ=0;
00084 double matchRatioPhi=0;
00085
00086 for (DTRecSegment4DCollection::const_iterator rechit = dtRecHits->begin(); rechit!=dtRecHits->end();++rechit) {
00087
00088 if ( !rechit->isValid()) continue;
00089
00090 double nhitsPhi = 0;
00091 double nhitsZ = 0;
00092
00093 if(rechit->hasZed()) {
00094 double countMuonDTHits = 0;
00095 double countAgreeingHits=0;
00096
00097 const DTRecSegment2D* segmZ;
00098 segmZ = dynamic_cast<const DTRecSegment2D*>(rechit->zSegment());
00099 nhitsZ = segmZ->recHits().size();
00100
00101 const vector<DTRecHit1D> hits1d = segmZ->specificRecHits();
00102 DTChamberId chamberSegIdT((segmZ->geographicalId()).rawId());
00103
00104
00105 for(trackingRecHit_iterator hit = dtHits.begin(); hit != dtHits.end(); ++hit) {
00106
00107 if ( !(*hit)->isValid()) continue;
00108
00109 DetId idT = (*hit)->geographicalId();
00110 DTChamberId dtDetIdHitT(idT.rawId());
00111 DTSuperLayerId dtDetLayerIdHitT(idT.rawId());
00112
00113 LocalPoint pointLocal = (*hit)->localPosition();
00114
00115 if ((chamberSegIdT==dtDetIdHitT) && (dtDetLayerIdHitT.superlayer()==2)) countMuonDTHits++;
00116
00117 for (vector<DTRecHit1D>::const_iterator hiti=hits1d.begin(); hiti!=hits1d.end(); hiti++) {
00118
00119 if ( !hiti->isValid()) continue;
00120
00121
00122 if(!(hiti->geographicalId().rawId()==idT.rawId())) continue;
00123
00124
00125 LocalPoint segLocal = hiti->localPosition();
00126
00127 if ((fabs(pointLocal.x()-segLocal.x())<ZCutParameter) &&
00128 (fabs(pointLocal.y()-segLocal.y())<ZCutParameter))
00129 countAgreeingHits++;
00130 }
00131 }
00132
00133 matchRatioZ = countMuonDTHits == 0 ? 0 : countAgreeingHits/countMuonDTHits;
00134 }
00135
00136 if(rechit->hasPhi()) {
00137 double countMuonDTHits = 0;
00138 double countAgreeingHits=0;
00139
00140
00141 const DTRecSegment2D* segmPhi;
00142 segmPhi = dynamic_cast<const DTRecSegment2D*>(rechit->phiSegment());
00143 nhitsPhi = segmPhi->recHits().size();
00144
00145 const vector<DTRecHit1D> hits1d = segmPhi->specificRecHits();
00146 DTChamberId chamberSegIdT((segmPhi->geographicalId()).rawId());
00147
00148
00149 for(trackingRecHit_iterator hit = dtHits.begin(); hit != dtHits.end(); ++hit) {
00150
00151 if ( !(*hit)->isValid()) continue;
00152
00153 DetId idT = (*hit)->geographicalId();
00154 DTChamberId dtDetIdHitT(idT.rawId());
00155 DTSuperLayerId dtDetLayerIdHitT(idT.rawId());
00156
00157 LocalPoint pointLocal = (*hit)->localPosition();
00158
00159 if ((chamberSegIdT==dtDetIdHitT)&&((dtDetLayerIdHitT.superlayer()==1)||(dtDetLayerIdHitT.superlayer()==3)))
00160 countMuonDTHits++;
00161
00162 for (vector<DTRecHit1D>::const_iterator hiti=hits1d.begin(); hiti!=hits1d.end(); hiti++) {
00163
00164 if ( !hiti->isValid()) continue;
00165
00166
00167 if(!(hiti->geographicalId().rawId()==idT.rawId())) continue;
00168
00169
00170 LocalPoint segLocal = hiti->localPosition();
00171
00172
00173 if ((fabs(pointLocal.x()-segLocal.x())<PhiCutParameter) &&
00174 (fabs(pointLocal.y()-segLocal.y())<PhiCutParameter))
00175 countAgreeingHits++;
00176 }
00177 }
00178
00179 matchRatioPhi = countMuonDTHits != 0 ? countAgreeingHits/countMuonDTHits : 0;
00180 }
00181
00182
00183
00184 if (dtTightMatch && nhitsPhi && nhitsZ) {
00185 if((matchRatioPhi>0.9)&&(matchRatioZ>0.9)) {
00186
00187 pointerTo4DSegments.push_back(&(*rechit));
00188 }
00189 } else {
00190 if((matchRatioPhi>0.9 && nhitsPhi)||(matchRatioZ>0.9 && nhitsZ)) {
00191
00192 pointerTo4DSegments.push_back(&(*rechit));
00193 }
00194 }
00195
00196 }
00197
00198 return pointerTo4DSegments;
00199 }
00200
00201
00202
00203 vector<const CSCSegment*> MuonSegmentMatcher::matchCSC(const reco::Track& muon, const edm::Event& event)
00204 {
00205
00206 using namespace edm;
00207
00208 edm::Handle<CSCSegmentCollection> allSegmentsCSC;
00209 event.getByLabel(CSCSegmentTags_, allSegmentsCSC);
00210
00211 vector<const CSCSegment*> pointerToCSCSegments;
00212
00213 double matchRatioCSC=0;
00214 int numCSC = 0;
00215 double CSCXCut = 0.001;
00216 double CSCYCut = 0.001;
00217 double countMuonCSCHits = 0;
00218
00219 for(CSCSegmentCollection::const_iterator segmentCSC = allSegmentsCSC->begin(); segmentCSC != allSegmentsCSC->end(); segmentCSC++) {
00220 double CSCcountAgreeingHits=0;
00221
00222 if ( !segmentCSC->isValid()) continue;
00223
00224 numCSC++;
00225 const vector<CSCRecHit2D>& CSCRechits2D = segmentCSC->specificRecHits();
00226 countMuonCSCHits = 0;
00227 CSCDetId myChamber((*segmentCSC).geographicalId().rawId());
00228
00229 for(trackingRecHit_iterator hitC = muon.recHitsBegin(); hitC != muon.recHitsEnd(); ++hitC) {
00230 if (!(*hitC)->isValid()) continue;
00231 if ( (*hitC)->geographicalId().det() != DetId::Muon ) continue;
00232 if ( (*hitC)->geographicalId().subdetId() != MuonSubdetId::CSC ) continue;
00233 if (!(*hitC)->isValid()) continue;
00234
00235
00236 DetId id = (*hitC)->geographicalId();
00237 CSCDetId cscDetIdHit(id.rawId());
00238
00239 if(!(cscDetIdHit.ring()==myChamber.ring())) continue;
00240 if(!(cscDetIdHit.station()==myChamber.station())) continue;
00241 if(!(cscDetIdHit.endcap()==myChamber.endcap())) continue;
00242 if(!(cscDetIdHit.chamber()==myChamber.chamber())) continue;
00243
00244 countMuonCSCHits++;
00245
00246 LocalPoint positionLocalCSC = (*hitC)->localPosition();
00247
00248 for (vector<CSCRecHit2D>::const_iterator hiti=CSCRechits2D.begin(); hiti!=CSCRechits2D.end(); hiti++) {
00249
00250 if ( !hiti->isValid()) continue;
00251 CSCDetId cscDetId((hiti->geographicalId()).rawId());
00252
00253 if ((*hitC)->geographicalId().rawId()!=(hiti->geographicalId()).rawId()) continue;
00254
00255 LocalPoint segLocalCSC = hiti->localPosition();
00256
00257
00258 if((fabs(positionLocalCSC.x()-segLocalCSC.x())<CSCXCut) &&
00259 (fabs(positionLocalCSC.y()-segLocalCSC.y())<CSCYCut)) {
00260 CSCcountAgreeingHits++;
00261
00262 }
00263 }
00264 }
00265
00266 matchRatioCSC = countMuonCSCHits == 0 ? 0 : CSCcountAgreeingHits/countMuonCSCHits;
00267
00268 if ((matchRatioCSC>0.9) && ((countMuonCSCHits>1) || !cscTightMatch)) pointerToCSCSegments.push_back(&(*segmentCSC));
00269
00270 }
00271
00272 return pointerToCSCSegments;
00273
00274 }
00275
00276
00277