Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-02-05 23:50:57

0001 #include "Alignment/CommonAlignmentProducer/interface/AlignmentCSCTrackSelector.h"
0002 
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 
0006 #include "FWCore/Framework/interface/Event.h"
0007 
0008 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
0009 #include "DataFormats/TrackReco/interface/Track.h"
0010 #include "DataFormats/DetId/interface/DetId.h"
0011 #include "DataFormats/MuonDetId/interface/CSCDetId.h"
0012 
0013 // constructor ----------------------------------------------------------------
0014 
0015 AlignmentCSCTrackSelector::AlignmentCSCTrackSelector(const edm::ParameterSet& cfg)
0016     : m_src(cfg.getParameter<edm::InputTag>("src")),
0017       m_stationA(cfg.getParameter<int>("stationA")),
0018       m_stationB(cfg.getParameter<int>("stationB")),
0019       m_minHitsDT(cfg.getParameter<int>("minHitsDT")),
0020       m_minHitsPerStation(cfg.getParameter<int>("minHitsPerStation")),
0021       m_maxHitsPerStation(cfg.getParameter<int>("maxHitsPerStation")) {}
0022 
0023 void AlignmentCSCTrackSelector::fillPSetDescription(edm::ParameterSetDescription& desc) {
0024   desc.add<int>("stationA", 0);
0025   desc.add<int>("stationB", 0);
0026   desc.add<int>("minHitsDT", 0);
0027   desc.add<int>("minHitsPerStation", 0);
0028   desc.add<int>("maxHitsPerStation", 0);
0029 }
0030 
0031 // destructor -----------------------------------------------------------------
0032 
0033 AlignmentCSCTrackSelector::~AlignmentCSCTrackSelector() {}
0034 
0035 // do selection ---------------------------------------------------------------
0036 
0037 AlignmentCSCTrackSelector::Tracks AlignmentCSCTrackSelector::select(const Tracks& tracks, const edm::Event& evt) const {
0038   Tracks result;
0039 
0040   for (auto const& track : tracks) {
0041     int hitsOnStationA = 0;
0042     int hitsOnStationB = 0;
0043 
0044     for (auto const& hit : track->recHits()) {
0045       DetId id = hit->geographicalId();
0046 
0047       if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::DT) {
0048         if (m_stationA == 0)
0049           hitsOnStationA++;
0050         if (m_stationB == 0)
0051           hitsOnStationB++;
0052       } else if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::CSC) {
0053         CSCDetId cscid(id.rawId());
0054         int station = (cscid.endcap() == 1 ? 1 : -1) * cscid.station();
0055 
0056         if (station == m_stationA)
0057           hitsOnStationA++;
0058         if (station == m_stationB)
0059           hitsOnStationB++;
0060 
0061       }  // end if CSC
0062     }  // end loop over hits
0063 
0064     bool stationAokay;
0065     if (m_stationA == 0)
0066       stationAokay = (m_minHitsDT <= hitsOnStationA);
0067     else
0068       stationAokay = (m_minHitsPerStation <= hitsOnStationA && hitsOnStationA <= m_maxHitsPerStation);
0069 
0070     bool stationBokay;
0071     if (m_stationB == 0)
0072       stationBokay = (m_minHitsDT <= hitsOnStationB);
0073     else
0074       stationBokay = (m_minHitsPerStation <= hitsOnStationB && hitsOnStationB <= m_maxHitsPerStation);
0075 
0076     if (stationAokay && stationBokay) {
0077       result.push_back(track);
0078     }
0079   }  // end loop over tracks
0080 
0081   return result;
0082 }