Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:14:07

0001 #include "L1Trigger/VertexFinder/interface/VertexProducer.h"
0002 #include "DataFormats/L1TrackTrigger/interface/TTTypes.h"
0003 #include "DataFormats/L1Trigger/interface/Vertex.h"
0004 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "Geometry/Records/interface/TrackerTopologyRcd.h"
0011 #include "L1Trigger/VertexFinder/interface/AlgoSettings.h"
0012 #include "L1Trigger/VertexFinder/interface/RecoVertex.h"
0013 #include "L1Trigger/VertexFinder/interface/VertexFinder.h"
0014 
0015 using namespace l1tVertexFinder;
0016 using namespace std;
0017 
0018 VertexProducer::VertexProducer(const edm::ParameterSet& iConfig)
0019     : l1TracksToken_(consumes<TTTrackCollectionView>(iConfig.getParameter<edm::InputTag>("l1TracksInputTag"))),
0020       tTopoToken(esConsumes<TrackerTopology, TrackerTopologyRcd>()),
0021       outputCollectionName_(iConfig.getParameter<std::string>("l1VertexCollectionName")),
0022       settings_(AlgoSettings(iConfig)) {
0023   // Get configuration parameters
0024 
0025   switch (settings_.vx_algo()) {
0026     case Algorithm::fastHisto:
0027       edm::LogInfo("VertexProducer") << "VertexProducer::Finding vertices using the fastHisto binning algorithm";
0028       break;
0029     case Algorithm::fastHistoEmulation:
0030       edm::LogInfo("VertexProducer")
0031           << "VertexProducer::Finding vertices using the emulation version of the fastHisto binning algorithm";
0032       break;
0033     case Algorithm::fastHistoLooseAssociation:
0034       edm::LogInfo("VertexProducer")
0035           << "VertexProducer::Finding vertices using the fastHistoLooseAssociation binning algorithm";
0036       break;
0037     case Algorithm::GapClustering:
0038       edm::LogInfo("VertexProducer") << "VertexProducer::Finding vertices using a gap clustering algorithm";
0039       break;
0040     case Algorithm::agglomerativeHierarchical:
0041       edm::LogInfo("VertexProducer") << "VertexProducer::Finding vertices using a Simple Merge Clustering algorithm";
0042       break;
0043     case Algorithm::DBSCAN:
0044       edm::LogInfo("VertexProducer") << "VertexProducer::Finding vertices using a DBSCAN algorithm";
0045       break;
0046     case Algorithm::PVR:
0047       edm::LogInfo("VertexProducer") << "VertexProducer::Finding vertices using a PVR algorithm";
0048       break;
0049     case Algorithm::adaptiveVertexReconstruction:
0050       edm::LogInfo("VertexProducer")
0051           << "VertexProducer::Finding vertices using an AdaptiveVertexReconstruction algorithm";
0052       break;
0053     case Algorithm::HPV:
0054       edm::LogInfo("VertexProducer") << "VertexProducer::Finding vertices using a Highest Pt Vertex algorithm";
0055       break;
0056     case Algorithm::Kmeans:
0057       edm::LogInfo("VertexProducer") << "VertexProducer::Finding vertices using a kmeans algorithm";
0058       break;
0059   }
0060 
0061   //--- Define EDM output to be written to file (if required)
0062   if (settings_.vx_algo() == Algorithm::fastHistoEmulation) {
0063     produces<l1t::VertexWordCollection>(outputCollectionName_ + "Emulation");
0064   } else {
0065     produces<l1t::VertexCollection>(outputCollectionName_);
0066   }
0067 }
0068 
0069 void VertexProducer::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0070   edm::Handle<TTTrackCollectionView> l1TracksHandle;
0071   iEvent.getByToken(l1TracksToken_, l1TracksHandle);
0072 
0073   std::vector<l1tVertexFinder::L1Track> l1Tracks;
0074   l1Tracks.reserve(l1TracksHandle->size());
0075   if (settings_.debug() > 1) {
0076     edm::LogInfo("VertexProducer") << "produce::Processing " << l1TracksHandle->size() << " tracks";
0077   }
0078   for (const auto& track : l1TracksHandle->ptrs()) {
0079     auto l1track = L1Track(track);
0080     // Check the minimum pT of the tracks
0081     // This is left here because it represents the smallest pT to be sent by the track finding boards
0082     // This has less to do with the algorithms than the constraints of what will be sent to the vertexing algorithm
0083     if (l1track.pt() >= settings_.vx_TrackMinPt()) {
0084       l1Tracks.push_back(l1track);
0085     } else {
0086       if (settings_.debug() > 2) {
0087         edm::LogInfo("VertexProducer") << "produce::Removing track with too low of a pt (" << l1track.pt() << ")\n"
0088                                        << "         word = " << l1track.getTTTrackPtr()->getTrackWord().to_string(2);
0089       }
0090     }
0091   }
0092   if (settings_.debug() > 1) {
0093     edm::LogInfo("VertexProducer") << "produce::Processing " << l1Tracks.size() << " tracks after minimum pt cut of"
0094                                    << settings_.vx_TrackMinPt() << " GeV";
0095   }
0096 
0097   VertexFinder vf(l1Tracks, settings_);
0098 
0099   switch (settings_.vx_algo()) {
0100     case Algorithm::fastHisto: {
0101       const TrackerTopology& tTopo = iSetup.getData(tTopoToken);
0102       vf.fastHisto(&tTopo);
0103       break;
0104     }
0105     case Algorithm::fastHistoEmulation:
0106       vf.fastHistoEmulation();
0107       break;
0108     case Algorithm::fastHistoLooseAssociation:
0109       vf.fastHistoLooseAssociation();
0110       break;
0111     case Algorithm::GapClustering:
0112       vf.GapClustering();
0113       break;
0114     case Algorithm::agglomerativeHierarchical:
0115       vf.agglomerativeHierarchicalClustering();
0116       break;
0117     case Algorithm::DBSCAN:
0118       vf.DBSCAN();
0119       break;
0120     case Algorithm::PVR:
0121       vf.PVR();
0122       break;
0123     case Algorithm::adaptiveVertexReconstruction:
0124       vf.adaptiveVertexReconstruction();
0125       break;
0126     case Algorithm::HPV:
0127       vf.HPV();
0128       break;
0129     case Algorithm::Kmeans:
0130       vf.Kmeans();
0131       break;
0132   }
0133 
0134   vf.sortVerticesInPt();
0135   vf.findPrimaryVertex();
0136 
0137   // //=== Store output EDM track and hardware stub collections.
0138   if (settings_.vx_algo() == Algorithm::fastHistoEmulation) {
0139     std::unique_ptr<l1t::VertexWordCollection> product_emulation =
0140         std::make_unique<l1t::VertexWordCollection>(vf.verticesEmulation().begin(), vf.verticesEmulation().end());
0141     iEvent.put(std::move(product_emulation), outputCollectionName_ + "Emulation");
0142   } else {
0143     std::unique_ptr<l1t::VertexCollection> product(new std::vector<l1t::Vertex>());
0144     for (const auto& vtx : vf.vertices()) {
0145       product->emplace_back(vtx.vertex());
0146     }
0147     iEvent.put(std::move(product), outputCollectionName_);
0148   }
0149 }
0150 DEFINE_FWK_MODULE(VertexProducer);