CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_6/src/RecoMuon/TrackingTools/src/SegmentsTrackAssociator.cc

Go to the documentation of this file.
00001 
00002 
00003 /*
00004  *  See header file for a description of this class.
00005  *
00006  *  $Date: 2009/06/26 15:18:25 $
00007  *  $Revision: 1.3 $
00008  *  \author C. Botta, G. Mila - INFN Torino
00009  */
00010 
00011 
00012 #include "RecoMuon/TrackingTools/interface/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   // The segment collections
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   //loop over recHit
00078   for(trackingRecHit_iterator recHit =  Track.recHitsBegin(); recHit != Track.recHitsEnd(); ++recHit){
00079 
00080     if(!(*recHit)->isValid()) continue;
00081 
00082     
00083     numRecHit++;
00084  
00085     //get the detector Id
00086     DetId idRivHit = (*recHit)->geographicalId();
00087       
00088 
00089     // DT recHits
00090     if (idRivHit.det() == DetId::Muon && idRivHit.subdetId() == MuonSubdetId::DT ) {
00091 
00092       // get the RecHit Local Position
00093       LocalPoint posTrackRecHit = (*recHit)->localPosition(); 
00094 
00095       DTRecSegment4DCollection::range range;    
00096       numRecHitDT++;
00097 
00098       // get the chamber Id
00099       DTChamberId chamberId(idRivHit.rawId());
00100       // get the segments of the chamber
00101       range = dtSegments->get(chamberId);
00102 
00103       // loop over segments
00104       for (segment = range.first; segment!=range.second; segment++){
00105 
00106         DetId id = segment->geographicalId();
00107         const GeomDet* det = theTrackingGeometry->idToDet(id);
00108         vector<const TrackingRecHit*> segments2D = (&(*segment))->recHits();
00109         // container for 4D segment recHits
00110         vector <const TrackingRecHit*> dtRecHits;
00111 
00112         for(vector<const TrackingRecHit*>::const_iterator segm2D = segments2D.begin();
00113           segm2D != segments2D.end();
00114           segm2D++) {
00115 
00116           vector <const TrackingRecHit*> rHits1D = (*segm2D)->recHits();
00117           for (int hit=0; hit<int(rHits1D.size()); hit++){
00118             dtRecHits.push_back(rHits1D[hit]);
00119           }
00120 
00121         }
00122         
00123         // loop over the recHit checking if there's the recHit of the track
00124         for (unsigned int hit = 0; hit < dtRecHits.size(); hit++) {       
00125           
00126           DetId idRivHitSeg = (*dtRecHits[hit]).geographicalId();
00127           LocalPoint posDTSegment=  segment->localPosition();
00128           LocalPoint posSegDTRecHit = (*dtRecHits[hit]).localPosition(); 
00129                     
00130           double rDT=sqrt(pow((posSegDTRecHit.x()-posTrackRecHit.x()),2) +pow((posSegDTRecHit.y()-posTrackRecHit.y()),2) + pow((posSegDTRecHit.z()-posTrackRecHit.z()),2));
00131           
00132           if (idRivHit == idRivHitSeg && rDT<0.0001){  
00133 
00134             if (selectedDtSegments.empty()){
00135                 selectedDtSegments.push_back(MuonTransientTrackingRecHit::specificBuild(det,&*segment));
00136                 LogTrace(metname) <<"First selected segment (from DT). Position : "<<posDTSegment<<"  Chamber : "<<segment->chamberId();
00137             }
00138             else{
00139               int check=0;
00140               for(int segm=0; segm < int(selectedDtSegments.size()); ++segm) {
00141                 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));
00142                 if(dist>30) check++;
00143               }     
00144                 
00145               if(check==int(selectedDtSegments.size())){
00146                 selectedDtSegments.push_back(MuonTransientTrackingRecHit::specificBuild(det,&*segment));
00147                 LogTrace(metname) <<"New DT selected segment. Position : "<<posDTSegment<<"  Chamber : "<<segment->chamberId();
00148               }      
00149             }   
00150           } // check to tag the segment as "selected"
00151 
00152         } // loop over segment recHits
00153 
00154       } // loop over DT segments
00155     }
00156     
00157    
00158     // CSC recHits
00159     if (idRivHit.det() == DetId::Muon && idRivHit.subdetId() == MuonSubdetId::CSC ) {
00160 
00161       // get the RecHit Local Position
00162       LocalPoint posTrackRecHit = (*recHit)->localPosition(); 
00163 
00164       CSCSegmentCollection::range range; 
00165       numRecHitCSC++;
00166 
00167       // get the chamber Id
00168       CSCDetId tempchamberId(idRivHit.rawId());
00169       
00170       int ring = tempchamberId.ring();
00171       int station = tempchamberId.station();
00172       int endcap = tempchamberId.endcap();
00173       int chamber = tempchamberId.chamber();    
00174       CSCDetId chamberId(endcap, station, ring, chamber, 0);
00175       
00176       // get the segments of the chamber
00177       range = cscSegments->get(chamberId);
00178       // loop over segments
00179       for(segment2 = range.first; segment2!=range.second; segment2++){
00180 
00181         DetId id2 = segment2->geographicalId();
00182         const GeomDet* det2 = theTrackingGeometry->idToDet(id2);
00183         
00184         // container for CSC segment recHits
00185         vector<const TrackingRecHit*> cscRecHits = (&(*segment2))->recHits();
00186         
00187         // loop over the recHit checking if there's the recHit of the track     
00188         for (unsigned int hit = 0; hit < cscRecHits.size(); hit++) { 
00189   
00190           DetId idRivHitSeg = (*cscRecHits[hit]).geographicalId();
00191           LocalPoint posSegCSCRecHit = (*cscRecHits[hit]).localPosition(); 
00192           LocalPoint posCSCSegment=  segment2->localPosition();
00193                   
00194           double rCSC=sqrt(pow((posSegCSCRecHit.x()-posTrackRecHit.x()),2) +pow((posSegCSCRecHit.y()-posTrackRecHit.y()),2) + pow((posSegCSCRecHit.z()-posTrackRecHit.z()),2));
00195 
00196           if (idRivHit==idRivHitSeg && rCSC < 0.0001){
00197 
00198             if (selectedCscSegments.empty()){
00199               selectedCscSegments.push_back(MuonTransientTrackingRecHit::specificBuild(det2,&*segment2));
00200               LogTrace(metname) <<"First selected segment (from CSC). Position: "<<posCSCSegment;
00201             }
00202             else{
00203               int check=0;
00204               for(int n=0; n< int(selectedCscSegments.size()); n++){
00205                 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));
00206                 if(dist>30) check++;
00207               }
00208               if(check==int(selectedCscSegments.size())){  
00209                 selectedCscSegments.push_back(MuonTransientTrackingRecHit::specificBuild(det2,&*segment2));
00210                 LogTrace(metname) <<"New CSC segment selected. Position : "<<posCSCSegment;     
00211               }
00212             }
00213             
00214           } // check to tag the segment as "selected"
00215             
00216         } // loop over segment recHits
00217         
00218       } // loop over DT segments
00219     } 
00220       
00221   } // loop over track recHits
00222     
00223   LogTrace(metname) <<"N recHit:"<<numRecHit;
00224   numRecHit=0;
00225   LogTrace(metname) <<"N recHit DT:"<<numRecHitDT;
00226   numRecHitDT=0;
00227   LogTrace(metname) <<"N recHit CSC:"<<numRecHitCSC;
00228   numRecHitCSC=0;
00229   
00230   copy(selectedDtSegments.begin(), selectedDtSegments.end(), back_inserter(selectedSegments));
00231   LogTrace(metname) <<"N selected Dt segments:"<<selectedDtSegments.size();
00232   copy(selectedCscSegments.begin(), selectedCscSegments.end(), back_inserter(selectedSegments));
00233   LogTrace(metname) <<"N selected Csc segments:"<<selectedCscSegments.size();
00234   LogTrace(metname) <<"N selected segments:"<<selectedSegments.size();
00235 
00236   return selectedSegments;
00237   
00238 }