![]() |
![]() |
00001 #include "Alignment/CommonAlignmentProducer/interface/AlignmentCSCTrackSelector.h" 00002 00003 #include "FWCore/MessageLogger/interface/MessageLogger.h" 00004 #include "FWCore/ParameterSet/interface/ParameterSet.h" 00005 00006 #include "FWCore/Framework/interface/Event.h" 00007 00008 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h" 00009 #include "DataFormats/TrackReco/interface/Track.h" 00010 #include "DataFormats/DetId/interface/DetId.h" 00011 #include "DataFormats/MuonDetId/interface/CSCDetId.h" 00012 00013 00014 // constructor ---------------------------------------------------------------- 00015 00016 AlignmentCSCTrackSelector::AlignmentCSCTrackSelector(const edm::ParameterSet & cfg) 00017 : m_src(cfg.getParameter<edm::InputTag>("src")) 00018 , m_stationA(cfg.getParameter<int>("stationA")) 00019 , m_stationB(cfg.getParameter<int>("stationB")) 00020 , m_minHitsDT(cfg.getParameter<int>("minHitsDT")) 00021 , m_minHitsPerStation(cfg.getParameter<int>("minHitsPerStation")) 00022 , m_maxHitsPerStation(cfg.getParameter<int>("maxHitsPerStation")) 00023 {} 00024 00025 // destructor ----------------------------------------------------------------- 00026 00027 AlignmentCSCTrackSelector::~AlignmentCSCTrackSelector() 00028 {} 00029 00030 00031 // do selection --------------------------------------------------------------- 00032 00033 AlignmentCSCTrackSelector::Tracks 00034 AlignmentCSCTrackSelector::select(const Tracks& tracks, const edm::Event& evt) const 00035 { 00036 Tracks result; 00037 00038 for (Tracks::const_iterator track = tracks.begin(); track != tracks.end(); ++track) { 00039 int hitsOnStationA = 0; 00040 int hitsOnStationB = 0; 00041 00042 for (trackingRecHit_iterator hit = (*track)->recHitsBegin(); hit != (*track)->recHitsEnd(); ++hit) { 00043 DetId id = (*hit)->geographicalId(); 00044 00045 if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::DT) { 00046 if (m_stationA == 0) hitsOnStationA++; 00047 if (m_stationB == 0) hitsOnStationB++; 00048 } 00049 else if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::CSC) { 00050 CSCDetId cscid(id.rawId()); 00051 int station = (cscid.endcap() == 1 ? 1 : -1) * cscid.station(); 00052 00053 if (station == m_stationA) hitsOnStationA++; 00054 if (station == m_stationB) hitsOnStationB++; 00055 00056 } // end if CSC 00057 } // end loop over hits 00058 00059 bool stationAokay; 00060 if (m_stationA == 0) stationAokay = (m_minHitsDT <= hitsOnStationA); 00061 else stationAokay = (m_minHitsPerStation <= hitsOnStationA && hitsOnStationA <= m_maxHitsPerStation); 00062 00063 bool stationBokay; 00064 if (m_stationB == 0) stationBokay = (m_minHitsDT <= hitsOnStationB); 00065 else stationBokay = (m_minHitsPerStation <= hitsOnStationB && hitsOnStationB <= m_maxHitsPerStation); 00066 00067 if (stationAokay && stationBokay) { 00068 result.push_back(*track); 00069 } 00070 } // end loop over tracks 00071 00072 return result; 00073 } 00074