File indexing completed on 2024-04-06 11:56:16
0001 #include "Alignment/CommonAlignmentProducer/interface/AlignmentCSCBeamHaloSelector.h"
0002
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005
0006 #include "DataFormats/DetId/interface/DetId.h"
0007 #include "DataFormats/MuonDetId/interface/MuonSubdetId.h"
0008 #include "DataFormats/MuonDetId/interface/CSCDetId.h"
0009 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
0010 #include "DataFormats/TrackingRecHit/interface/TrackingRecHitFwd.h"
0011
0012
0013
0014 AlignmentCSCBeamHaloSelector::AlignmentCSCBeamHaloSelector(const edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)
0015 : m_minStations(iConfig.getParameter<unsigned int>("minStations")),
0016 m_minHitsPerStation(iConfig.getParameter<unsigned int>("minHitsPerStation")) {
0017 edm::LogInfo("AlignmentCSCBeamHaloSelector")
0018 << "Acceptable tracks must have at least " << m_minHitsPerStation << " hits in " << m_minStations
0019 << " different CSC stations." << std::endl;
0020 }
0021
0022
0023
0024 AlignmentCSCBeamHaloSelector::~AlignmentCSCBeamHaloSelector() {}
0025
0026
0027
0028 AlignmentCSCBeamHaloSelector::Tracks AlignmentCSCBeamHaloSelector::select(const Tracks &tracks,
0029 const edm::Event &iEvent) const {
0030 Tracks result;
0031
0032 for (auto const &track : tracks) {
0033 std::map<int, unsigned int> station_map;
0034
0035 for (auto const &hit : track->recHits()) {
0036 DetId id = hit->geographicalId();
0037 if (id.det() == DetId::Muon && id.subdetId() == MuonSubdetId::CSC) {
0038 CSCDetId cscid(id.rawId());
0039 int station = (cscid.endcap() == 1 ? 1 : -1) * cscid.station();
0040
0041 std::map<int, unsigned int>::const_iterator station_iter = station_map.find(station);
0042 if (station_iter == station_map.end()) {
0043 station_map[station] = 0;
0044 }
0045 station_map[station]++;
0046 }
0047 }
0048
0049 unsigned int stations = 0;
0050 for (std::map<int, unsigned int>::const_iterator station_iter = station_map.begin();
0051 station_iter != station_map.end();
0052 ++station_iter) {
0053 if (station_iter->second > m_minHitsPerStation)
0054 stations++;
0055 }
0056 if (stations >= m_minStations) {
0057 result.push_back(track);
0058 }
0059 }
0060
0061 return result;
0062 }