Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:59:48

0001 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
0002 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
0003 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
0004 
0005 #include "CalibTracker/SiStripCommon/interface/SiStripOnTrackClusterTableProducerBase.h"
0006 
0007 SiStripOnTrackClusterTableProducerBase::~SiStripOnTrackClusterTableProducerBase() {}
0008 
0009 namespace {
0010   int findTrackIndex(const edm::View<reco::Track>& tracks, const reco::Track* track) {
0011     for (auto iTr = tracks.begin(); iTr != tracks.end(); ++iTr) {
0012       if (&(*iTr) == track) {
0013         return iTr - tracks.begin();
0014       }
0015     }
0016     return -2;
0017   }
0018 }  // namespace
0019 
0020 void SiStripOnTrackClusterTableProducerBase::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0021   edm::Handle<edm::View<reco::Track>> tracks;
0022   iEvent.getByToken(m_tracks_token, tracks);
0023   edm::Handle<TrajTrackAssociationCollection> trajTrackAssociations;
0024   iEvent.getByToken(m_association_token, trajTrackAssociations);
0025 
0026   std::vector<OnTrackCluster> clusters{};
0027 
0028   for (const auto& assoc : *trajTrackAssociations) {
0029     const auto traj = assoc.key.get();
0030     const auto track = assoc.val.get();
0031 
0032     for (const auto& meas : traj->measurements()) {
0033       const auto& trajState = meas.updatedState();
0034       if (!trajState.isValid())
0035         continue;
0036 
0037       // there can be 2 (stereo module), 1 (no stereo module), or 0 (no strip hit) clusters per measurement
0038       const auto trechit = meas.recHit()->hit();
0039       const auto simple1d = dynamic_cast<const SiStripRecHit1D*>(trechit);
0040       const auto simple = dynamic_cast<const SiStripRecHit2D*>(trechit);
0041       const auto matched = dynamic_cast<const SiStripMatchedRecHit2D*>(trechit);
0042       if (matched) {
0043         clusters.emplace_back(matched->monoId(), &matched->monoCluster(), traj, track, meas);
0044         clusters.emplace_back(matched->stereoId(), &matched->stereoCluster(), traj, track, meas);
0045       } else if (simple) {
0046         clusters.emplace_back(simple->geographicalId().rawId(), simple->cluster().get(), traj, track, meas);
0047       } else if (simple1d) {
0048         clusters.emplace_back(simple1d->geographicalId().rawId(), simple1d->cluster().get(), traj, track, meas);
0049       }
0050     }
0051   }
0052 
0053   auto out = std::make_unique<nanoaod::FlatTable>(clusters.size(), m_name, false, m_extension);
0054   if (!m_extension) {
0055     std::vector<int> c_trackindex;
0056     c_trackindex.reserve(clusters.size());
0057     std::vector<uint32_t> c_rawid;
0058     c_rawid.reserve(clusters.size());
0059     for (const auto clus : clusters) {
0060       c_trackindex.push_back(findTrackIndex(*tracks, clus.track));
0061       c_rawid.push_back(clus.det);
0062     }
0063     addColumn(out.get(), "trackindex", c_trackindex, "Track index");
0064     addColumn(out.get(), "rawid", c_rawid, "DetId");
0065   }
0066   fillTable(clusters, *tracks, out.get(), iSetup);
0067   iEvent.put(std::move(out));
0068 }