00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 #include "DQMOffline/Muon/src/SegmentsTrackAssociator.h"
00013
00014 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00015 #include "DataFormats/TrackReco/interface/Track.h"
00016 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00017 #include "DataFormats/MuonDetId/interface/MuonSubdetId.h"
00018 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00019 #include "DataFormats/DetId/interface/DetId.h"
00020 #include "DataFormats/Common/interface/getRef.h"
00021 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00022 #include "DataFormats/DTRecHit/interface/DTRecSegment4D.h"
00023 #include "DataFormats/TrackingRecHit/interface/RecSegment.h"
00024
00025 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
00026 #include "DataFormats/CSCRecHit/interface/CSCSegmentCollection.h"
00027 #include "DataFormats/RPCRecHit/interface/RPCRecHitCollection.h"
00028 #include "DataFormats/DTRecHit/interface/DTRecHitCollection.h"
00029 #include "DataFormats/CSCRecHit/interface/CSCRecHit2DCollection.h"
00030
00031 #include "MagneticField/Engine/interface/MagneticField.h"
00032 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00033 #include "Geometry/Records/interface/GlobalTrackingGeometryRecord.h"
00034 #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
00035 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00036
00037 #include <vector>
00038
00039 using namespace edm;
00040 using namespace std;
00041
00042
00043
00044 SegmentsTrackAssociator::SegmentsTrackAssociator(const ParameterSet& iConfig)
00045 {
00046 theDTSegmentLabel = iConfig.getUntrackedParameter<InputTag>("segmentsDt");
00047 theCSCSegmentLabel = iConfig.getUntrackedParameter<InputTag>("segmentsCSC");
00048 theSegmentContainerName = iConfig.getUntrackedParameter<InputTag>("SelectedSegments");
00049
00050 numRecHit=0;
00051 numRecHitDT=0;
00052 numRecHitCSC=0;
00053 metname = "SegmentsTrackAssociator";
00054 }
00055
00056
00057 SegmentsTrackAssociator::~SegmentsTrackAssociator();
00058
00059
00060 MuonTransientTrackingRecHit::MuonRecHitContainer SegmentsTrackAssociator::associate(const Event& iEvent, const EventSetup& iSetup, const reco::Track& Track){
00061
00062
00063 Handle<DTRecSegment4DCollection> dtSegments;
00064 iEvent.getByLabel(theDTSegmentLabel, dtSegments);
00065 Handle<CSCSegmentCollection> cscSegments;
00066 iEvent.getByLabel(theCSCSegmentLabel, cscSegments);
00067 ESHandle<GlobalTrackingGeometry> theTrackingGeometry;
00068 iSetup.get<GlobalTrackingGeometryRecord>().get(theTrackingGeometry);
00069
00070
00071 DTRecSegment4DCollection::const_iterator segment;
00072 CSCSegmentCollection::const_iterator segment2;
00073 MuonTransientTrackingRecHit::MuonRecHitContainer selectedSegments;
00074 MuonTransientTrackingRecHit::MuonRecHitContainer selectedDtSegments;
00075 MuonTransientTrackingRecHit::MuonRecHitContainer selectedCscSegments;
00076
00077
00078 for(trackingRecHit_iterator recHit = Track.recHitsBegin(); recHit != Track.recHitsEnd(); ++recHit){
00079
00080 numRecHit++;
00081
00082
00083 DetId idRivHit = (*recHit)->geographicalId();
00084 LocalPoint posTrackRecHit = (*recHit)->localPosition();
00085
00086
00087
00088 if (idRivHit.det() == DetId::Muon && idRivHit.subdetId() == MuonSubdetId::DT ) {
00089
00090 DTRecSegment4DCollection::range range;
00091 numRecHitDT++;
00092
00093 DTChamberId chamberId(idRivHit.rawId());
00094
00095 range = dtSegments->get(chamberId);
00096
00097
00098 for (segment = range.first; segment!=range.second; segment++){
00099
00100 DetId id = segment->geographicalId();
00101 const GeomDet* det = theTrackingGeometry->idToDet(id);
00102 vector<const TrackingRecHit*> segments2D = (&(*segment))->recHits();
00103
00104 vector <const TrackingRecHit*> dtRecHits;
00105
00106 for(vector<const TrackingRecHit*>::const_iterator segm2D = segments2D.begin();
00107 segm2D != segments2D.end();
00108 segm2D++) {
00109
00110 vector <const TrackingRecHit*> rHits1D = (*segm2D)->recHits();
00111 for (int hit=0; hit<int(rHits1D.size()); hit++){
00112 dtRecHits.push_back(rHits1D[hit]);
00113 }
00114
00115 }
00116
00117
00118 for (unsigned int hit = 0; hit < dtRecHits.size(); hit++) {
00119
00120 DetId idRivHitSeg = (*dtRecHits[hit]).geographicalId();
00121 LocalPoint posDTSegment= segment->localPosition();
00122 LocalPoint posSegDTRecHit = (*dtRecHits[hit]).localPosition();
00123
00124 double rDT=sqrt(pow((posSegDTRecHit.x()-posTrackRecHit.x()),2) +pow((posSegDTRecHit.y()-posTrackRecHit.y()),2) + pow((posSegDTRecHit.z()-posTrackRecHit.z()),2));
00125
00126 if (idRivHit == idRivHitSeg && rDT<0.0001){
00127
00128 if (selectedDtSegments.empty()){
00129 selectedDtSegments.push_back(MuonTransientTrackingRecHit::specificBuild(det,&*segment));
00130 LogTrace(metname) <<"First selected segment (from DT). Position : "<<posDTSegment<<" Chamber : "<<segment->chamberId();
00131 }
00132 else{
00133 int check=0;
00134 for(int segm=0; segm < int(selectedDtSegments.size()); ++segm) {
00135 double dist=sqrt(pow((((*(selectedDtSegments[segm])).localPosition()).x()-posDTSegment.x()),2) +pow((((*(selectedDtSegments[segm])).localPosition()).y()-posDTSegment.y()),2) + pow((((*(selectedDtSegments[segm])).localPosition()).z()-posDTSegment.z()),2));
00136 if(dist>30) check++;
00137 }
00138
00139 if(check==int(selectedDtSegments.size())){
00140 selectedDtSegments.push_back(MuonTransientTrackingRecHit::specificBuild(det,&*segment));
00141 LogTrace(metname) <<"New DT selected segment. Position : "<<posDTSegment<<" Chamber : "<<segment->chamberId();
00142 }
00143 }
00144 }
00145
00146 }
00147
00148 }
00149 }
00150
00151
00152
00153 if (idRivHit.det() == DetId::Muon && idRivHit.subdetId() == MuonSubdetId::CSC ) {
00154
00155 CSCSegmentCollection::range range;
00156 numRecHitCSC++;
00157
00158
00159 CSCDetId tempchamberId(idRivHit.rawId());
00160
00161 int ring = tempchamberId.ring();
00162 int station = tempchamberId.station();
00163 int endcap = tempchamberId.endcap();
00164 int chamber = tempchamberId.chamber();
00165 CSCDetId chamberId(endcap, station, ring, chamber, 0);
00166
00167
00168 range = cscSegments->get(chamberId);
00169
00170 for(segment2 = range.first; segment2!=range.second; segment2++){
00171
00172 DetId id2 = segment2->geographicalId();
00173 const GeomDet* det2 = theTrackingGeometry->idToDet(id2);
00174
00175
00176 vector<const TrackingRecHit*> cscRecHits = (&(*segment2))->recHits();
00177
00178
00179 for (unsigned int hit = 0; hit < cscRecHits.size(); hit++) {
00180
00181 DetId idRivHitSeg = (*cscRecHits[hit]).geographicalId();
00182 LocalPoint posSegCSCRecHit = (*cscRecHits[hit]).localPosition();
00183 LocalPoint posCSCSegment= segment2->localPosition();
00184
00185 double rCSC=sqrt(pow((posSegCSCRecHit.x()-posTrackRecHit.x()),2) +pow((posSegCSCRecHit.y()-posTrackRecHit.y()),2) + pow((posSegCSCRecHit.z()-posTrackRecHit.z()),2));
00186
00187 if (idRivHit==idRivHitSeg && rCSC < 0.0001){
00188
00189 if (selectedCscSegments.empty()){
00190 selectedCscSegments.push_back(MuonTransientTrackingRecHit::specificBuild(det2,&*segment2));
00191 LogTrace(metname) <<"First selected segment (from CSC). Position: "<<posCSCSegment;
00192 }
00193 else{
00194 int check=0;
00195 for(int n=0; n< int(selectedCscSegments.size()); n++){
00196 double dist = sqrt(pow(((*(selectedCscSegments[n])).localPosition().x()-posCSCSegment.x()),2) +pow(((*(selectedCscSegments[n])).localPosition().y()-posCSCSegment.y()),2) + pow(((*(selectedCscSegments[n])).localPosition().z()-posCSCSegment.z()),2));
00197 if(dist>30) check++;
00198 }
00199 if(check==int(selectedCscSegments.size())){
00200 selectedCscSegments.push_back(MuonTransientTrackingRecHit::specificBuild(det2,&*segment2));
00201 LogTrace(metname) <<"New CSC segment selected. Position : "<<posCSCSegment;
00202 }
00203 }
00204
00205 }
00206
00207 }
00208
00209 }
00210 }
00211
00212 }
00213
00214 LogTrace(metname) <<"N recHit:"<<numRecHit;
00215 numRecHit=0;
00216 LogTrace(metname) <<"N recHit DT:"<<numRecHitDT;
00217 numRecHitDT=0;
00218 LogTrace(metname) <<"N recHit CSC:"<<numRecHitCSC;
00219 numRecHitCSC=0;
00220
00221 copy(selectedDtSegments.begin(), selectedDtSegments.end(), back_inserter(selectedSegments));
00222 LogTrace(metname) <<"N selected Dt segments:"<<selectedDtSegments.size();
00223 copy(selectedCscSegments.begin(), selectedCscSegments.end(), back_inserter(selectedSegments));
00224 LogTrace(metname) <<"N selected Csc segments:"<<selectedCscSegments.size();
00225 LogTrace(metname) <<"N selected segments:"<<selectedSegments.size();
00226
00227 return selectedSegments;
00228
00229 }