Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:27:19

0001 
0002 
0003 /*
0004  *  See header file for a description of this class.
0005  *
0006  *  \author C. Botta, G. Mila - INFN Torino
0007  */
0008 
0009 #include "RecoMuon/TrackingTools/interface/SegmentsTrackAssociator.h"
0010 
0011 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
0012 #include "DataFormats/TrackReco/interface/Track.h"
0013 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
0014 #include "DataFormats/MuonDetId/interface/MuonSubdetId.h"
0015 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
0016 #include "DataFormats/DetId/interface/DetId.h"
0017 #include "DataFormats/Common/interface/getRef.h"
0018 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0019 #include "DataFormats/DTRecHit/interface/DTRecSegment4D.h"
0020 #include "DataFormats/TrackingRecHit/interface/RecSegment.h"
0021 
0022 #include "DataFormats/RPCRecHit/interface/RPCRecHitCollection.h"
0023 #include "DataFormats/DTRecHit/interface/DTRecHitCollection.h"
0024 #include "DataFormats/CSCRecHit/interface/CSCRecHit2DCollection.h"
0025 
0026 #include "MagneticField/Engine/interface/MagneticField.h"
0027 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0028 #include "Geometry/Records/interface/GlobalTrackingGeometryRecord.h"
0029 #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
0030 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0031 
0032 #include <vector>
0033 
0034 using namespace edm;
0035 using namespace std;
0036 
0037 SegmentsTrackAssociator::SegmentsTrackAssociator(const ParameterSet& iConfig, edm::ConsumesCollector& iC) {
0038   theDTSegmentLabel = iConfig.getUntrackedParameter<InputTag>("segmentsDt");
0039   theCSCSegmentLabel = iConfig.getUntrackedParameter<InputTag>("segmentsCSC");
0040   theSegmentContainerName = iConfig.getUntrackedParameter<InputTag>("SelectedSegments");
0041 
0042   numRecHit = 0;
0043   numRecHitDT = 0;
0044   numRecHitCSC = 0;
0045   metname = "SegmentsTrackAssociator";
0046 
0047   dtSegmentsToken = iC.consumes<DTRecSegment4DCollection>(theDTSegmentLabel);
0048   cscSegmentsToken = iC.consumes<CSCSegmentCollection>(theCSCSegmentLabel);
0049   trackingGeometryToken = iC.esConsumes();
0050 }
0051 
0052 SegmentsTrackAssociator::~SegmentsTrackAssociator() {}
0053 
0054 MuonTransientTrackingRecHit::MuonRecHitContainer SegmentsTrackAssociator::associate(const Event& iEvent,
0055                                                                                     const EventSetup& iSetup,
0056                                                                                     const reco::Track& Track) {
0057   // The segment collections
0058   Handle<DTRecSegment4DCollection> dtSegments;
0059   iEvent.getByToken(dtSegmentsToken, dtSegments);
0060   Handle<CSCSegmentCollection> cscSegments;
0061   iEvent.getByToken(cscSegmentsToken, cscSegments);
0062   ESHandle<GlobalTrackingGeometry> theTrackingGeometry = iSetup.getHandle(trackingGeometryToken);
0063 
0064   DTRecSegment4DCollection::const_iterator segment;
0065   CSCSegmentCollection::const_iterator segment2;
0066   MuonTransientTrackingRecHit::MuonRecHitContainer selectedSegments;
0067   MuonTransientTrackingRecHit::MuonRecHitContainer selectedDtSegments;
0068   MuonTransientTrackingRecHit::MuonRecHitContainer selectedCscSegments;
0069 
0070   //loop over recHit
0071   for (auto const& recHit : Track.recHits()) {
0072     if (!recHit->isValid())
0073       continue;
0074 
0075     numRecHit++;
0076 
0077     //get the detector Id
0078     DetId idRivHit = recHit->geographicalId();
0079 
0080     // DT recHits
0081     if (idRivHit.det() == DetId::Muon && idRivHit.subdetId() == MuonSubdetId::DT) {
0082       // get the RecHit Local Position
0083       LocalPoint posTrackRecHit = recHit->localPosition();
0084 
0085       DTRecSegment4DCollection::range range;
0086       numRecHitDT++;
0087 
0088       // get the chamber Id
0089       DTChamberId chamberId(idRivHit.rawId());
0090       // get the segments of the chamber
0091       range = dtSegments->get(chamberId);
0092 
0093       // loop over segments
0094       for (segment = range.first; segment != range.second; segment++) {
0095         DetId id = segment->geographicalId();
0096         const GeomDet* det = theTrackingGeometry->idToDet(id);
0097         vector<const TrackingRecHit*> segments2D = (&(*segment))->recHits();
0098         // container for 4D segment recHits
0099         vector<const TrackingRecHit*> dtRecHits;
0100 
0101         for (vector<const TrackingRecHit*>::const_iterator segm2D = segments2D.begin(); segm2D != segments2D.end();
0102              segm2D++) {
0103           vector<const TrackingRecHit*> rHits1D = (*segm2D)->recHits();
0104           for (int hit = 0; hit < int(rHits1D.size()); hit++) {
0105             dtRecHits.push_back(rHits1D[hit]);
0106           }
0107         }
0108 
0109         // loop over the recHit checking if there's the recHit of the track
0110         for (unsigned int hit = 0; hit < dtRecHits.size(); hit++) {
0111           DetId idRivHitSeg = (*dtRecHits[hit]).geographicalId();
0112           LocalPoint posDTSegment = segment->localPosition();
0113           LocalPoint posSegDTRecHit = (*dtRecHits[hit]).localPosition();
0114 
0115           double rDT = sqrt(pow((posSegDTRecHit.x() - posTrackRecHit.x()), 2) +
0116                             pow((posSegDTRecHit.y() - posTrackRecHit.y()), 2) +
0117                             pow((posSegDTRecHit.z() - posTrackRecHit.z()), 2));
0118 
0119           if (idRivHit == idRivHitSeg && rDT < 0.0001) {
0120             if (selectedDtSegments.empty()) {
0121               selectedDtSegments.push_back(MuonTransientTrackingRecHit::specificBuild(det, &*segment));
0122               LogTrace(metname) << "First selected segment (from DT). Position : " << posDTSegment
0123                                 << "  Chamber : " << segment->chamberId();
0124             } else {
0125               int check = 0;
0126               for (int segm = 0; segm < int(selectedDtSegments.size()); ++segm) {
0127                 double dist = sqrt(pow((((*(selectedDtSegments[segm])).localPosition()).x() - posDTSegment.x()), 2) +
0128                                    pow((((*(selectedDtSegments[segm])).localPosition()).y() - posDTSegment.y()), 2) +
0129                                    pow((((*(selectedDtSegments[segm])).localPosition()).z() - posDTSegment.z()), 2));
0130                 if (dist > 30)
0131                   check++;
0132               }
0133 
0134               if (check == int(selectedDtSegments.size())) {
0135                 selectedDtSegments.push_back(MuonTransientTrackingRecHit::specificBuild(det, &*segment));
0136                 LogTrace(metname) << "New DT selected segment. Position : " << posDTSegment
0137                                   << "  Chamber : " << segment->chamberId();
0138               }
0139             }
0140           }  // check to tag the segment as "selected"
0141 
0142         }  // loop over segment recHits
0143 
0144       }  // loop over DT segments
0145     }
0146 
0147     // CSC recHits
0148     if (idRivHit.det() == DetId::Muon && idRivHit.subdetId() == MuonSubdetId::CSC) {
0149       // get the RecHit Local Position
0150       LocalPoint posTrackRecHit = recHit->localPosition();
0151 
0152       CSCSegmentCollection::range range;
0153       numRecHitCSC++;
0154 
0155       // get the chamber Id
0156       CSCDetId tempchamberId(idRivHit.rawId());
0157 
0158       int ring = tempchamberId.ring();
0159       int station = tempchamberId.station();
0160       int endcap = tempchamberId.endcap();
0161       int chamber = tempchamberId.chamber();
0162       CSCDetId chamberId(endcap, station, ring, chamber, 0);
0163 
0164       // get the segments of the chamber
0165       range = cscSegments->get(chamberId);
0166       // loop over segments
0167       for (segment2 = range.first; segment2 != range.second; segment2++) {
0168         DetId id2 = segment2->geographicalId();
0169         const GeomDet* det2 = theTrackingGeometry->idToDet(id2);
0170 
0171         // container for CSC segment recHits
0172         vector<const TrackingRecHit*> cscRecHits = (&(*segment2))->recHits();
0173 
0174         // loop over the recHit checking if there's the recHit of the track
0175         for (unsigned int hit = 0; hit < cscRecHits.size(); hit++) {
0176           DetId idRivHitSeg = (*cscRecHits[hit]).geographicalId();
0177           LocalPoint posSegCSCRecHit = (*cscRecHits[hit]).localPosition();
0178           LocalPoint posCSCSegment = segment2->localPosition();
0179 
0180           double rCSC = sqrt(pow((posSegCSCRecHit.x() - posTrackRecHit.x()), 2) +
0181                              pow((posSegCSCRecHit.y() - posTrackRecHit.y()), 2) +
0182                              pow((posSegCSCRecHit.z() - posTrackRecHit.z()), 2));
0183 
0184           if (idRivHit == idRivHitSeg && rCSC < 0.0001) {
0185             if (selectedCscSegments.empty()) {
0186               selectedCscSegments.push_back(MuonTransientTrackingRecHit::specificBuild(det2, &*segment2));
0187               LogTrace(metname) << "First selected segment (from CSC). Position: " << posCSCSegment;
0188             } else {
0189               int check = 0;
0190               for (int n = 0; n < int(selectedCscSegments.size()); n++) {
0191                 double dist = sqrt(pow(((*(selectedCscSegments[n])).localPosition().x() - posCSCSegment.x()), 2) +
0192                                    pow(((*(selectedCscSegments[n])).localPosition().y() - posCSCSegment.y()), 2) +
0193                                    pow(((*(selectedCscSegments[n])).localPosition().z() - posCSCSegment.z()), 2));
0194                 if (dist > 30)
0195                   check++;
0196               }
0197               if (check == int(selectedCscSegments.size())) {
0198                 selectedCscSegments.push_back(MuonTransientTrackingRecHit::specificBuild(det2, &*segment2));
0199                 LogTrace(metname) << "New CSC segment selected. Position : " << posCSCSegment;
0200               }
0201             }
0202 
0203           }  // check to tag the segment as "selected"
0204 
0205         }  // loop over segment recHits
0206 
0207       }  // loop over DT segments
0208     }
0209 
0210   }  // loop over track recHits
0211 
0212   LogTrace(metname) << "N recHit:" << numRecHit;
0213   numRecHit = 0;
0214   LogTrace(metname) << "N recHit DT:" << numRecHitDT;
0215   numRecHitDT = 0;
0216   LogTrace(metname) << "N recHit CSC:" << numRecHitCSC;
0217   numRecHitCSC = 0;
0218 
0219   copy(selectedDtSegments.begin(), selectedDtSegments.end(), back_inserter(selectedSegments));
0220   LogTrace(metname) << "N selected Dt segments:" << selectedDtSegments.size();
0221   copy(selectedCscSegments.begin(), selectedCscSegments.end(), back_inserter(selectedSegments));
0222   LogTrace(metname) << "N selected Csc segments:" << selectedCscSegments.size();
0223   LogTrace(metname) << "N selected segments:" << selectedSegments.size();
0224 
0225   return selectedSegments;
0226 }