Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-06-04 04:35:10

0001 //////////////////////////
0002 //  Producer by Anders  //
0003 //     and Emmanuele    //
0004 //    july 2012 @ CU    //
0005 //////////////////////////
0006 
0007 ////////////////////
0008 // FRAMEWORK HEADERS
0009 #include "FWCore/PluginManager/interface/ModuleDef.h"
0010 #include "FWCore/Framework/interface/MakerMacros.h"
0011 //
0012 #include "FWCore/Framework/interface/one/EDProducer.h"
0013 #include "FWCore/Framework/interface/Event.h"
0014 #include "FWCore/Framework/interface/ESHandle.h"
0015 #include "FWCore/Framework/interface/EventSetup.h"
0016 #include "FWCore/ParameterSet/interface/FileInPath.h"
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018 //
0019 #include "FWCore/Utilities/interface/InputTag.h"
0020 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0021 #include "FWCore/ServiceRegistry/interface/Service.h"
0022 
0023 ///////////////////////
0024 // DATA FORMATS HEADERS
0025 #include "DataFormats/Common/interface/Handle.h"
0026 #include "DataFormats/Common/interface/Ref.h"
0027 //
0028 #include "DataFormats/DetId/interface/DetId.h"
0029 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
0030 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
0031 #include "DataFormats/Common/interface/DetSetVector.h"
0032 #include "DataFormats/L1TrackTrigger/interface/TTBV.h"
0033 //
0034 #include "L1Trigger/TrackFindingTracklet/interface/SLHCEvent.h"
0035 #include "L1Trigger/TrackFindingTracklet/interface/L1TStub.h"
0036 
0037 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
0038 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
0039 #include "SimDataFormats/Track/interface/SimTrack.h"
0040 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
0041 #include "SimDataFormats/Vertex/interface/SimVertex.h"
0042 #include "SimDataFormats/Vertex/interface/SimVertexContainer.h"
0043 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
0044 #include "SimTracker/TrackTriggerAssociation/interface/TTClusterAssociationMap.h"
0045 //
0046 #include "DataFormats/Math/interface/LorentzVector.h"
0047 #include "DataFormats/Math/interface/Vector3D.h"
0048 //
0049 #include "DataFormats/L1TrackTrigger/interface/TTCluster.h"
0050 #include "DataFormats/L1TrackTrigger/interface/TTStub.h"
0051 #include "DataFormats/L1TrackTrigger/interface/TTTrack.h"
0052 #include "DataFormats/L1TrackTrigger/interface/TTTrack_TrackWord.h"
0053 #include "DataFormats/L1TrackTrigger/interface/TTTypes.h"
0054 #include "DataFormats/L1TrackTrigger/interface/TTDTC.h"
0055 #include "L1Trigger/TrackTrigger/interface/Setup.h"
0056 //
0057 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0058 #include "DataFormats/Candidate/interface/Candidate.h"
0059 //
0060 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementPoint.h"
0061 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0062 #include "DataFormats/SiPixelDetId/interface/PixelChannelIdentifier.h"
0063 #include "TrackingTools/GeomPropagators/interface/HelixArbitraryPlaneCrossing.h"
0064 
0065 ////////////////////////////
0066 // DETECTOR GEOMETRY HEADERS
0067 #include "MagneticField/Engine/interface/MagneticField.h"
0068 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0069 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0070 #include "Geometry/CommonTopologies/interface/PixelGeomDetUnit.h"
0071 #include "Geometry/CommonTopologies/interface/PixelGeomDetType.h"
0072 #include "Geometry/TrackerGeometryBuilder/interface/PixelTopologyBuilder.h"
0073 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0074 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
0075 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0076 //
0077 #include "Geometry/Records/interface/StackedTrackerGeometryRecord.h"
0078 
0079 ///////////////
0080 // Tracklet emulation
0081 #include "L1Trigger/TrackFindingTracklet/interface/Settings.h"
0082 #include "L1Trigger/TrackFindingTracklet/interface/Sector.h"
0083 #include "L1Trigger/TrackFindingTracklet/interface/Track.h"
0084 #include "L1Trigger/TrackFindingTracklet/interface/TrackletEventProcessor.h"
0085 #include "L1Trigger/TrackFindingTracklet/interface/ChannelAssignment.h"
0086 #include "L1Trigger/TrackFindingTracklet/interface/Tracklet.h"
0087 #include "L1Trigger/TrackFindingTracklet/interface/Residual.h"
0088 #include "L1Trigger/TrackFindingTracklet/interface/Stub.h"
0089 #include "L1Trigger/TrackFindingTracklet/interface/StubKiller.h"
0090 #include "L1Trigger/TrackFindingTracklet/interface/StubStreamData.h"
0091 #include "L1Trigger/TrackFindingTracklet/interface/HitPatternHelper.h"
0092 
0093 ////////////////
0094 // PHYSICS TOOLS
0095 #include "RecoTracker/TkSeedGenerator/interface/FastHelix.h"
0096 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
0097 
0098 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementVector.h"
0099 #include "DataFormats/GeometrySurface/interface/BoundPlane.h"
0100 
0101 #include "L1Trigger/TrackTrigger/interface/StubPtConsistency.h"
0102 #include "L1Trigger/TrackTrigger/interface/L1TrackQuality.h"
0103 
0104 //////////////
0105 // STD HEADERS
0106 #include <memory>
0107 #include <string>
0108 #include <iostream>
0109 #include <fstream>
0110 
0111 //////////////
0112 // NAMESPACES
0113 using namespace edm;
0114 using namespace std;
0115 using namespace tt;
0116 using namespace trklet;
0117 
0118 //////////////////////////////
0119 //                          //
0120 //     CLASS DEFINITION     //
0121 //                          //
0122 //////////////////////////////
0123 
0124 /////////////////////////////////////
0125 // this class is needed to make a map
0126 // between different types of stubs
0127 struct L1TStubCompare {
0128 public:
0129   bool operator()(const trklet::L1TStub& a, const trklet::L1TStub& b) const {
0130     if (a.x() != b.x())
0131       return (b.x() > a.x());
0132     else if (a.y() != b.y())
0133       return (b.y() > a.y());
0134     else if (a.z() != b.z())
0135       return (a.z() > b.z());
0136     else
0137       return a.bend() > b.bend();
0138   }
0139 };
0140 
0141 class L1FPGATrackProducer : public edm::one::EDProducer<edm::one::WatchRuns> {
0142 public:
0143   /// Constructor/destructor
0144   explicit L1FPGATrackProducer(const edm::ParameterSet& iConfig);
0145   ~L1FPGATrackProducer() override;
0146 
0147 private:
0148   int eventnum;
0149 
0150   /// Containers of parameters passed by python configuration file
0151   edm::ParameterSet config;
0152 
0153   bool readMoreMcTruth_;
0154 
0155   /// File path for configuration files
0156   edm::FileInPath fitPatternFile;
0157   edm::FileInPath memoryModulesFile;
0158   edm::FileInPath processingModulesFile;
0159   edm::FileInPath wiresFile;
0160 
0161   edm::FileInPath tableTEDFile;
0162   edm::FileInPath tableTREFile;
0163 
0164   string asciiEventOutName_;
0165   std::ofstream asciiEventOut_;
0166 
0167   // settings containing various constants for the tracklet processing
0168   trklet::Settings settings_;
0169 
0170   // event processor for the tracklet track finding
0171   trklet::TrackletEventProcessor eventProcessor;
0172 
0173   // used to "kill" stubs from a selected area of the detector
0174   StubKiller* stubKiller_;
0175   int failScenario_;
0176 
0177   unsigned int nHelixPar_;
0178   bool extended_;
0179   bool reduced_;
0180 
0181   bool trackQuality_;
0182   std::unique_ptr<L1TrackQuality> trackQualityModel_;
0183 
0184   std::map<string, vector<int>> dtclayerdisk;
0185 
0186   edm::InputTag MCTruthClusterInputTag;
0187   edm::InputTag MCTruthStubInputTag;
0188   edm::InputTag TrackingParticleInputTag;
0189 
0190   const edm::EDGetTokenT<reco::BeamSpot> getTokenBS_;
0191   const edm::EDGetTokenT<TTDTC> getTokenDTC_;
0192   edm::EDGetTokenT<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>> getTokenTTClusterMCTruth_;
0193   edm::EDGetTokenT<std::vector<TrackingParticle>> getTokenTrackingParticle_;
0194 
0195   // ED output token for clock and bit accurate tracks
0196   const edm::EDPutTokenT<Streams> putTokenTracks_;
0197   // ED output token for clock and bit accurate stubs
0198   const edm::EDPutTokenT<StreamsStub> putTokenStubs_;
0199   // ChannelAssignment token
0200   const ESGetToken<ChannelAssignment, ChannelAssignmentRcd> esGetTokenChannelAssignment_;
0201   // helper class to assign tracks to channel
0202   const ChannelAssignment* channelAssignment_;
0203 
0204   // helper class to store DTC configuration
0205   const Setup* setup_;
0206   // helper class to store configuration needed by HitPatternHelper
0207   const hph::Setup* setupHPH_;
0208 
0209   // Setup token
0210   const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> esGetTokenBfield_;
0211   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> esGetTokenTGeom_;
0212   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> esGetTokenTTopo_;
0213   const edm::ESGetToken<tt::Setup, tt::SetupRcd> esGetToken_;
0214   const edm::ESGetToken<hph::Setup, hph::SetupRcd> esGetTokenHPH_;
0215 
0216   /// ///////////////// ///
0217   /// MANDATORY METHODS ///
0218   void beginRun(const edm::Run& run, const edm::EventSetup& iSetup) override;
0219   void endRun(edm::Run const&, edm::EventSetup const&) override;
0220   void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override;
0221 };
0222 
0223 //////////////
0224 // CONSTRUCTOR
0225 L1FPGATrackProducer::L1FPGATrackProducer(edm::ParameterSet const& iConfig)
0226     : config(iConfig),
0227       readMoreMcTruth_(iConfig.getParameter<bool>("readMoreMcTruth")),
0228       MCTruthClusterInputTag(readMoreMcTruth_ ? config.getParameter<edm::InputTag>("MCTruthClusterInputTag")
0229                                               : edm::InputTag()),
0230       MCTruthStubInputTag(readMoreMcTruth_ ? config.getParameter<edm::InputTag>("MCTruthStubInputTag")
0231                                            : edm::InputTag()),
0232       TrackingParticleInputTag(readMoreMcTruth_ ? iConfig.getParameter<edm::InputTag>("TrackingParticleInputTag")
0233                                                 : edm::InputTag()),
0234       // book ED products
0235       getTokenBS_(consumes<reco::BeamSpot>(config.getParameter<edm::InputTag>("BeamSpotSource"))),
0236       getTokenDTC_(consumes<TTDTC>(edm::InputTag(iConfig.getParameter<edm::InputTag>("InputTagTTDTC")))),
0237       // book ED output token for clock and bit accurate tracks
0238       putTokenTracks_(produces<Streams>("Level1TTTracks")),
0239       // book ED output token for clock and bit accurate stubs
0240       putTokenStubs_(produces<StreamsStub>("Level1TTTracks")),
0241       // book ES products
0242       esGetTokenChannelAssignment_(esConsumes<ChannelAssignment, ChannelAssignmentRcd, Transition::BeginRun>()),
0243       esGetTokenBfield_(esConsumes<edm::Transition::BeginRun>()),
0244       esGetTokenTGeom_(esConsumes()),
0245       esGetTokenTTopo_(esConsumes()),
0246       esGetToken_(esConsumes<tt::Setup, tt::SetupRcd, edm::Transition::BeginRun>()),
0247       esGetTokenHPH_(esConsumes<hph::Setup, hph::SetupRcd, edm::Transition::BeginRun>()) {
0248   if (readMoreMcTruth_) {
0249     getTokenTTClusterMCTruth_ = consumes<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>>(MCTruthClusterInputTag);
0250     getTokenTrackingParticle_ = consumes<std::vector<TrackingParticle>>(TrackingParticleInputTag);
0251   }
0252 
0253   produces<std::vector<TTTrack<Ref_Phase2TrackerDigi_>>>("Level1TTTracks").setBranchAlias("Level1TTTracks");
0254 
0255   asciiEventOutName_ = iConfig.getUntrackedParameter<string>("asciiFileName", "");
0256 
0257   fitPatternFile = iConfig.getParameter<edm::FileInPath>("fitPatternFile");
0258   processingModulesFile = iConfig.getParameter<edm::FileInPath>("processingModulesFile");
0259   memoryModulesFile = iConfig.getParameter<edm::FileInPath>("memoryModulesFile");
0260   wiresFile = iConfig.getParameter<edm::FileInPath>("wiresFile");
0261 
0262   failScenario_ = iConfig.getUntrackedParameter<int>("FailScenario", 0);
0263 
0264   extended_ = iConfig.getParameter<bool>("Extended");
0265   reduced_ = iConfig.getParameter<bool>("Reduced");
0266   nHelixPar_ = iConfig.getParameter<unsigned int>("Hnpar");
0267 
0268   if (extended_) {
0269     tableTEDFile = iConfig.getParameter<edm::FileInPath>("tableTEDFile");
0270     tableTREFile = iConfig.getParameter<edm::FileInPath>("tableTREFile");
0271   }
0272 
0273   // initial ES products
0274   channelAssignment_ = nullptr;
0275   setup_ = nullptr;
0276 
0277   // --------------------------------------------------------------------------------
0278   // set options in Settings based on inputs from configuration files
0279   // --------------------------------------------------------------------------------
0280 
0281   settings_.setExtended(extended_);
0282   settings_.setReduced(reduced_);
0283   settings_.setNHelixPar(nHelixPar_);
0284 
0285   settings_.setFitPatternFile(fitPatternFile.fullPath());
0286   settings_.setProcessingModulesFile(processingModulesFile.fullPath());
0287   settings_.setMemoryModulesFile(memoryModulesFile.fullPath());
0288   settings_.setWiresFile(wiresFile.fullPath());
0289 
0290   settings_.setFakefit(iConfig.getParameter<bool>("Fakefit"));
0291   settings_.setStoreTrackBuilderOutput(iConfig.getParameter<bool>("StoreTrackBuilderOutput"));
0292   settings_.setRemovalType(iConfig.getParameter<string>("RemovalType"));
0293   settings_.setDoMultipleMatches(iConfig.getParameter<bool>("DoMultipleMatches"));
0294 
0295   if (extended_) {
0296     settings_.setTableTEDFile(tableTEDFile.fullPath());
0297     settings_.setTableTREFile(tableTREFile.fullPath());
0298 
0299     //FIXME: The TED and TRE tables are currently disabled by default, so we
0300     //need to allow for the additional tracklets that will eventually be
0301     //removed by these tables, once they are finalized
0302     settings_.setNbitstrackletindex(15);
0303   }
0304 
0305   eventnum = 0;
0306   if (not asciiEventOutName_.empty()) {
0307     asciiEventOut_.open(asciiEventOutName_.c_str());
0308   }
0309 
0310   if (settings_.debugTracklet()) {
0311     edm::LogVerbatim("Tracklet") << "fit pattern :     " << fitPatternFile.fullPath()
0312                                  << "\n process modules : " << processingModulesFile.fullPath()
0313                                  << "\n memory modules :  " << memoryModulesFile.fullPath()
0314                                  << "\n wires          :  " << wiresFile.fullPath();
0315     if (extended_) {
0316       edm::LogVerbatim("Tracklet") << "table_TED    :  " << tableTEDFile.fullPath()
0317                                    << "\n table_TRE    :  " << tableTREFile.fullPath();
0318     }
0319   }
0320 
0321   trackQuality_ = iConfig.getParameter<bool>("TrackQuality");
0322   if (trackQuality_) {
0323     trackQualityModel_ = std::make_unique<L1TrackQuality>(iConfig.getParameter<edm::ParameterSet>("TrackQualityPSet"));
0324   }
0325   if (settings_.storeTrackBuilderOutput() && (settings_.doMultipleMatches() || !settings_.removalType().empty())) {
0326     cms::Exception exception("ConfigurationNotSupported.");
0327     exception.addContext("L1FPGATrackProducer::produce");
0328     if (settings_.doMultipleMatches())
0329       exception << "Storing of TrackBuilder output does not support doMultipleMatches.";
0330     if (!settings_.removalType().empty())
0331       exception << "Storing of TrackBuilder output does not support duplicate removal.";
0332     throw exception;
0333   }
0334 }
0335 
0336 /////////////
0337 // DESTRUCTOR
0338 L1FPGATrackProducer::~L1FPGATrackProducer() {
0339   if (asciiEventOut_.is_open()) {
0340     asciiEventOut_.close();
0341   }
0342 }
0343 
0344 ///////END RUN
0345 //
0346 void L1FPGATrackProducer::endRun(const edm::Run& run, const edm::EventSetup& iSetup) {}
0347 
0348 ////////////
0349 // BEGIN JOB
0350 void L1FPGATrackProducer::beginRun(const edm::Run& run, const edm::EventSetup& iSetup) {
0351   ////////////////////////
0352   // GET MAGNETIC FIELD //
0353   const MagneticField* theMagneticField = &iSetup.getData(esGetTokenBfield_);
0354   double mMagneticFieldStrength = theMagneticField->inTesla(GlobalPoint(0, 0, 0)).z();
0355   settings_.setBfield(mMagneticFieldStrength);
0356 
0357   setup_ = &iSetup.getData(esGetToken_);
0358 
0359   settings_.passSetup(setup_);
0360 
0361   setupHPH_ = &iSetup.getData(esGetTokenHPH_);
0362   // Tracklet pattern reco output channel info.
0363   channelAssignment_ = &iSetup.getData(esGetTokenChannelAssignment_);
0364   // initialize the tracklet event processing (this sets all the processing & memory modules, wiring, etc)
0365   eventProcessor.init(settings_, setup_);
0366 }
0367 
0368 //////////
0369 // PRODUCE
0370 void L1FPGATrackProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0371   typedef std::map<trklet::L1TStub,
0372                    edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>,
0373                    L1TStubCompare>
0374       stubMapType;
0375   typedef std::map<unsigned int,
0376                    edm::Ref<edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_>>, TTStub<Ref_Phase2TrackerDigi_>>>
0377       stubIndexMapType;
0378   typedef edm::Ref<edmNew::DetSetVector<TTCluster<Ref_Phase2TrackerDigi_>>, TTCluster<Ref_Phase2TrackerDigi_>>
0379       TTClusterRef;
0380 
0381   /// Prepare output
0382   auto L1TkTracksForOutput = std::make_unique<std::vector<TTTrack<Ref_Phase2TrackerDigi_>>>();
0383 
0384   stubMapType stubMap;
0385   stubIndexMapType stubIndexMap;
0386 
0387   ////////////
0388   // GET BS //
0389   edm::Handle<reco::BeamSpot> beamSpotHandle;
0390   iEvent.getByToken(getTokenBS_, beamSpotHandle);
0391   math::XYZPoint bsPosition = beamSpotHandle->position();
0392 
0393   eventnum++;
0394   trklet::SLHCEvent ev;
0395   ev.setEventNum(eventnum);
0396   ev.setIP(bsPosition.x(), bsPosition.y());
0397 
0398   // tracking particles
0399   edm::Handle<std::vector<TrackingParticle>> TrackingParticleHandle;
0400   if (readMoreMcTruth_)
0401     iEvent.getByToken(getTokenTrackingParticle_, TrackingParticleHandle);
0402 
0403   // tracker topology
0404   const TrackerTopology* const tTopo = &iSetup.getData(esGetTokenTTopo_);
0405   const TrackerGeometry* const theTrackerGeom = &iSetup.getData(esGetTokenTGeom_);
0406 
0407   // check killing stubs for detector degradation studies
0408   // if failType = 0, StubKiller does not kill any modules
0409   int failType = 0;
0410   if (failScenario_ < 0 || failScenario_ > 9) {
0411     edm::LogVerbatim("Tracklet") << "Invalid fail scenario! Ignoring input";
0412   } else
0413     failType = failScenario_;
0414 
0415   stubKiller_ = new StubKiller();
0416   stubKiller_->initialise(failType, tTopo, theTrackerGeom);
0417 
0418   ////////////////////////
0419   // GET THE PRIMITIVES //
0420   edm::Handle<TTDTC> handleDTC;
0421   iEvent.getByToken<TTDTC>(getTokenDTC_, handleDTC);
0422 
0423   // must be defined for code to compile, even if it's not used unless readMoreMcTruth_ is true
0424   map<edm::Ptr<TrackingParticle>, int> translateTP;
0425 
0426   // MC truth association maps
0427   edm::Handle<TTClusterAssociationMap<Ref_Phase2TrackerDigi_>> MCTruthTTClusterHandle;
0428   if (readMoreMcTruth_) {
0429     iEvent.getByToken(getTokenTTClusterMCTruth_, MCTruthTTClusterHandle);
0430 
0431     ////////////////////////////////////////////////
0432     /// LOOP OVER TRACKING PARTICLES & GET SIMTRACKS
0433 
0434     int ntps = 1;  //count from 1 ; 0 will mean invalid
0435 
0436     int this_tp = 0;
0437     if (readMoreMcTruth_) {
0438       for (const auto& iterTP : *TrackingParticleHandle) {
0439         edm::Ptr<TrackingParticle> tp_ptr(TrackingParticleHandle, this_tp);
0440         this_tp++;
0441 
0442         // only keep TPs producing a cluster
0443         if (MCTruthTTClusterHandle->findTTClusterRefs(tp_ptr).empty())
0444           continue;
0445 
0446         if (iterTP.g4Tracks().empty()) {
0447           continue;
0448         }
0449 
0450         int sim_eventid = iterTP.g4Tracks().at(0).eventId().event();
0451         int sim_type = iterTP.pdgId();
0452         float sim_pt = iterTP.pt();
0453         float sim_eta = iterTP.eta();
0454         float sim_phi = iterTP.phi();
0455 
0456         float vx = iterTP.vertex().x();
0457         float vy = iterTP.vertex().y();
0458         float vz = iterTP.vertex().z();
0459 
0460         if (sim_pt < 1.0 || std::abs(vz) > 100.0 || hypot(vx, vy) > 50.0)
0461           continue;
0462 
0463         ev.addL1SimTrack(sim_eventid, ntps, sim_type, sim_pt, sim_eta, sim_phi, vx, vy, vz);
0464 
0465         translateTP[tp_ptr] = ntps;
0466         ntps++;
0467 
0468       }  //end loop over TPs
0469     }
0470 
0471   }  // end if (readMoreMcTruth_)
0472 
0473   /////////////////////////////////
0474   /// READ DTC STUB INFORMATION ///
0475   /////////////////////////////////
0476 
0477   // Process stubs in each region and channel within that tracking region
0478   unsigned int theStubIndex = 0;
0479   for (const int& region : handleDTC->tfpRegions()) {
0480     for (const int& channel : handleDTC->tfpChannels()) {
0481       // Get the DTC name & ID from the channel
0482       unsigned int atcaSlot = channel % 12;
0483       string dtcname = settings_.slotToDTCname(atcaSlot);
0484       if (channel % 24 >= 12)
0485         dtcname = "neg" + dtcname;
0486       dtcname += (channel < 24) ? "_A" : "_B";  // which detector region
0487       int dtcId = setup_->dtcId(region, channel);
0488 
0489       // Get the stubs from the DTC
0490       const tt::StreamStub& streamFromDTC{handleDTC->stream(region, channel)};
0491 
0492       // Prepare the DTC stubs for the IR
0493       for (size_t stubIndex = 0; stubIndex < streamFromDTC.size(); ++stubIndex) {
0494         const tt::FrameStub& stub{streamFromDTC[stubIndex]};
0495         const TTStubRef& stubRef = stub.first;
0496 
0497         if (stubRef.isNull())
0498           continue;
0499 
0500         const GlobalPoint& ttPos = setup_->stubPos(stubRef);
0501 
0502         //Get the 2 bits for the layercode
0503         string layerword = stub.second.to_string().substr(61, 2);
0504         unsigned int layercode = 2 * (layerword[0] - '0') + layerword[1] - '0';
0505         assert(layercode < 4);
0506 
0507         //translation from the two bit layercode to the layer/disk number of each of the
0508         //12 channels (dtcs)
0509         // FIX: take this from DTC cabling map.
0510         static const int layerdisktab[12][4] = {{0, 6, 8, 10},
0511                                                 {0, 7, 9, -1},
0512                                                 {1, 7, -1, -1},
0513                                                 {6, 8, 10, -1},
0514                                                 {2, 7, -1, -1},
0515                                                 {2, 9, -1, -1},
0516                                                 {3, 4, -1, -1},
0517                                                 {4, -1, -1, -1},
0518                                                 {5, -1, -1, -1},
0519                                                 {5, 8, -1, -1},
0520                                                 {6, 9, -1, -1},
0521                                                 {7, 10, -1, -1}};
0522 
0523         int layerdisk = layerdisktab[channel % 12][layercode];
0524         assert(layerdisk != -1);
0525 
0526         //Get the 36 bit word - skip the lowest 3 buts (status and layer code)
0527         constexpr int DTCLinkWordSize = 64;
0528         constexpr int StubWordSize = 36;
0529         constexpr int LayerandStatusCodeSize = 3;
0530         string stubword =
0531             stub.second.to_string().substr(DTCLinkWordSize - StubWordSize - LayerandStatusCodeSize, StubWordSize);
0532         string stubwordhex = "";
0533 
0534         //Loop over the 9 words in the 36 bit stub word
0535         for (unsigned int i = 0; i < 9; i++) {
0536           bitset<4> bits(stubword.substr(i * 4, 4));
0537           ulong val = bits.to_ulong();
0538           stubwordhex += ((val < 10) ? ('0' + val) : ('A' + val - 10));
0539         }
0540 
0541         /// Get the Inner and Outer TTCluster
0542         edm::Ref<edmNew::DetSetVector<TTCluster<Ref_Phase2TrackerDigi_>>, TTCluster<Ref_Phase2TrackerDigi_>>
0543             innerCluster = stub.first->clusterRef(0);
0544         edm::Ref<edmNew::DetSetVector<TTCluster<Ref_Phase2TrackerDigi_>>, TTCluster<Ref_Phase2TrackerDigi_>>
0545             outerCluster = stub.first->clusterRef(1);
0546 
0547         // -----------------------------------------------------
0548         // check module orientation, if flipped, need to store that information for track fit
0549         // -----------------------------------------------------
0550 
0551         const DetId innerDetId = innerCluster->getDetId();
0552         const GeomDetUnit* det_inner = theTrackerGeom->idToDetUnit(innerDetId);
0553         const auto* theGeomDet_inner = dynamic_cast<const PixelGeomDetUnit*>(det_inner);
0554         const PixelTopology* topol_inner = dynamic_cast<const PixelTopology*>(&(theGeomDet_inner->specificTopology()));
0555 
0556         MeasurementPoint coords_inner = innerCluster->findAverageLocalCoordinatesCentered();
0557         LocalPoint clustlp_inner = topol_inner->localPosition(coords_inner);
0558         GlobalPoint posStub_inner = theGeomDet_inner->surface().toGlobal(clustlp_inner);
0559 
0560         const DetId outerDetId = outerCluster->getDetId();
0561         const GeomDetUnit* det_outer = theTrackerGeom->idToDetUnit(outerDetId);
0562         const auto* theGeomDet_outer = dynamic_cast<const PixelGeomDetUnit*>(det_outer);
0563         const PixelTopology* topol_outer = dynamic_cast<const PixelTopology*>(&(theGeomDet_outer->specificTopology()));
0564 
0565         MeasurementPoint coords_outer = outerCluster->findAverageLocalCoordinatesCentered();
0566         LocalPoint clustlp_outer = topol_outer->localPosition(coords_outer);
0567         GlobalPoint posStub_outer = theGeomDet_outer->surface().toGlobal(clustlp_outer);
0568 
0569         bool isFlipped = (posStub_outer.mag() < posStub_inner.mag());
0570 
0571         vector<int> assocTPs;
0572 
0573         for (unsigned int iClus = 0; iClus <= 1; iClus++) {  // Loop over both clusters that make up stub.
0574 
0575           const TTClusterRef& ttClusterRef = stubRef->clusterRef(iClus);
0576 
0577           // Now identify all TP's contributing to either cluster in stub.
0578           if (readMoreMcTruth_) {
0579             vector<edm::Ptr<TrackingParticle>> vecTpPtr =
0580                 MCTruthTTClusterHandle->findTrackingParticlePtrs(ttClusterRef);
0581 
0582             for (const edm::Ptr<TrackingParticle>& tpPtr : vecTpPtr) {
0583               if (translateTP.find(tpPtr) != translateTP.end()) {
0584                 if (iClus == 0) {
0585                   assocTPs.push_back(translateTP.at(tpPtr));
0586                 } else {
0587                   assocTPs.push_back(-translateTP.at(tpPtr));
0588                 }
0589                 // N.B. Since not all tracking particles are stored in InputData::vTPs_, sometimes no match will be found.
0590               } else {
0591                 assocTPs.push_back(0);
0592               }
0593             }
0594           }
0595         }
0596 
0597         double stubbend = stubRef->bendFE();  //stubRef->rawBend()
0598         if (ttPos.z() < -120) {
0599           stubbend = -stubbend;
0600         }
0601 
0602         bool barrel = (layerdisk < N_LAYER);
0603         // See  https://github.com/cms-sw/cmssw/tree/master/Geometry/TrackerNumberingBuilder
0604         enum TypeBarrel { nonBarrel = 0, tiltedMinus = 1, tiltedPlus = 2, flat = 3 };
0605         const TypeBarrel type = static_cast<TypeBarrel>(tTopo->tobSide(innerDetId));
0606         bool tiltedBarrel = barrel && (type == tiltedMinus || type == tiltedPlus);
0607         unsigned int tiltedRingId = 0;
0608         // Tilted module ring no. (Increasing 1 to 12 as |z| increases).
0609         if (tiltedBarrel) {
0610           tiltedRingId = tTopo->tobRod(innerDetId);
0611           if (type == tiltedMinus) {
0612             unsigned int layp1 = 1 + layerdisk;  // Setup counts from 1
0613             unsigned int nTilted = setup_->numTiltedLayerRing(layp1);
0614             tiltedRingId = 1 + nTilted - tiltedRingId;
0615           }
0616         }
0617         // Endcap module ring number (1-15) in endcap disks.
0618         unsigned int endcapRingId = barrel ? 0 : tTopo->tidRing(innerDetId);
0619 
0620         const unsigned int intDetId = innerDetId.rawId();
0621 
0622         // check killing stubs for detector degredation studies
0623         const TTStub<Ref_Phase2TrackerDigi_>* theStub = &(*stubRef);
0624         bool killThisStub = stubKiller_->killStub(theStub);
0625         if (!killThisStub) {
0626           ev.addStub(dtcname,
0627                      region,
0628                      layerdisk,
0629                      stubwordhex,
0630                      setup_->psModule(dtcId),
0631                      isFlipped,
0632                      tiltedBarrel,
0633                      tiltedRingId,
0634                      endcapRingId,
0635                      intDetId,
0636                      ttPos.x(),
0637                      ttPos.y(),
0638                      ttPos.z(),
0639                      stubbend,
0640                      stubRef->innerClusterPosition(),
0641                      assocTPs,
0642                      theStubIndex);
0643 
0644           const trklet::L1TStub& lastStub = ev.lastStub();
0645           stubMap[lastStub] = stubRef;
0646           stubIndexMap[lastStub.uniqueIndex()] = stub.first;
0647           theStubIndex++;
0648         }
0649       }
0650     }
0651   }
0652 
0653   //////////////////////////
0654   // NOW RUN THE L1 tracking
0655 
0656   if (!asciiEventOutName_.empty()) {
0657     ev.write(asciiEventOut_);
0658   }
0659 
0660   const std::vector<trklet::Track>& tracks = eventProcessor.tracks();
0661 
0662   // max number of projection layers
0663   const unsigned int maxNumProjectionLayers = channelAssignment_->maxNumProjectionLayers();
0664   // number of track channels
0665   const unsigned int numStreamsTrack = N_SECTOR * channelAssignment_->numChannelsTrack();
0666   // number of stub channels
0667   const unsigned int numStreamsStub = N_SECTOR * channelAssignment_->numChannelsStub();
0668   // number of seeding layers
0669   const unsigned int numSeedingLayers = channelAssignment_->numSeedingLayers();
0670   // max number of stub channel per track
0671   const unsigned int numStubChannel = maxNumProjectionLayers + numSeedingLayers;
0672   // number of stub channels if all seed types streams padded to have same number of stub channels (for coding simplicity)
0673   const unsigned int numStreamsStubRaw = numStreamsTrack * numStubChannel;
0674 
0675   // Streams formatted to allow this code to run outside CMSSW.
0676   vector<vector<string>> streamsTrackRaw(numStreamsTrack);
0677   vector<vector<StubStreamData>> streamsStubRaw(numStreamsStubRaw);
0678 
0679   // this performs the actual tracklet event processing
0680   eventProcessor.event(ev, streamsTrackRaw, streamsStubRaw);
0681 
0682   for (const auto& track : tracks) {
0683     if (track.duplicate())
0684       continue;
0685 
0686     // this is where we create the TTTrack object
0687     double tmp_rinv = track.rinv(settings_);
0688     double tmp_phi = track.phi0(settings_);
0689     double tmp_tanL = track.tanL(settings_);
0690     double tmp_z0 = track.z0(settings_);
0691     double tmp_d0 = track.d0(settings_);
0692     double tmp_chi2rphi = track.chisqrphi();
0693     double tmp_chi2rz = track.chisqrz();
0694     unsigned int tmp_hit = track.hitpattern();
0695 
0696     TTTrack<Ref_Phase2TrackerDigi_> aTrack(tmp_rinv,
0697                                            tmp_phi,
0698                                            tmp_tanL,
0699                                            tmp_z0,
0700                                            tmp_d0,
0701                                            tmp_chi2rphi,
0702                                            tmp_chi2rz,
0703                                            0,
0704                                            0,
0705                                            0,
0706                                            tmp_hit,
0707                                            settings_.nHelixPar(),
0708                                            settings_.bfield());
0709 
0710     unsigned int trksector = track.sector();
0711     unsigned int trkseed = (unsigned int)abs(track.seed());
0712 
0713     aTrack.setPhiSector(trksector);
0714     aTrack.setTrackSeedType(trkseed);
0715 
0716     const vector<trklet::L1TStub>& stubptrs = track.stubs();
0717     vector<trklet::L1TStub> stubs;
0718 
0719     stubs.reserve(stubptrs.size());
0720     for (const auto& stubptr : stubptrs) {
0721       stubs.push_back(stubptr);
0722     }
0723 
0724     int countStubs = 0;
0725     stubMapType::const_iterator it;
0726     stubIndexMapType::const_iterator itIndex;
0727     for (const auto& itstubs : stubs) {
0728       itIndex = stubIndexMap.find(itstubs.uniqueIndex());
0729       if (itIndex != stubIndexMap.end()) {
0730         aTrack.addStubRef(itIndex->second);
0731         countStubs = countStubs + 1;
0732       } else {
0733         // could not find stub in stub map
0734       }
0735     }
0736 
0737     // pt consistency
0738     aTrack.setStubPtConsistency(
0739         StubPtConsistency::getConsistency(aTrack, theTrackerGeom, tTopo, settings_.bfield(), settings_.nHelixPar()));
0740 
0741     // set track word before TQ MVA calculated which uses track word variables
0742     aTrack.setTrackWordBits();
0743 
0744     if (trackQuality_) {
0745       trackQualityModel_->setL1TrackQuality(aTrack);
0746     }
0747 
0748     //    hph::HitPatternHelper hph(setupHPH_, tmp_hit, tmp_tanL, tmp_z0);
0749     //    if (trackQuality_) {
0750     //      trackQualityModel_->setBonusFeatures(hph.bonusFeatures());
0751     //    }
0752 
0753     // test track word
0754     //aTrack.testTrackWordBits();
0755 
0756     // set track word again to set MVA variable from TTTrack into track word
0757     aTrack.setTrackWordBits();
0758     // test track word
0759     //aTrack.testTrackWordBits();
0760 
0761     L1TkTracksForOutput->push_back(aTrack);
0762   }
0763 
0764   iEvent.put(std::move(L1TkTracksForOutput), "Level1TTTracks");
0765 
0766   // produce clock and bit accurate stream output tracks and stubs.
0767   // from end of tracklet pattern recognition.
0768   // Convertion here is from stream format that allows this code to run
0769   // outside CMSSW to the EDProduct one.
0770   Streams streamsTrack(numStreamsTrack);
0771   StreamsStub streamsStub(numStreamsStub);
0772 
0773   for (unsigned int chanTrk = 0; chanTrk < numStreamsTrack; chanTrk++) {
0774     for (unsigned int itk = 0; itk < streamsTrackRaw[chanTrk].size(); itk++) {
0775       std::string bitsTrk = streamsTrackRaw[chanTrk][itk];
0776       int iSeed = chanTrk % channelAssignment_->numChannelsTrack();  // seed type
0777       streamsTrack[chanTrk].emplace_back(bitsTrk);
0778 
0779       const unsigned int chanStubOffsetIn = chanTrk * numStubChannel;
0780       const unsigned int chanStubOffsetOut = channelAssignment_->offsetStub(chanTrk);
0781       const unsigned int numProjLayers = channelAssignment_->numProjectionLayers(iSeed);
0782       TTBV hitMap(0, numProjLayers + numSeedingLayers);
0783       // remove padding from stub stream
0784       for (unsigned int iproj = 0; iproj < numStubChannel; iproj++) {
0785         // FW current has one (perhaps invalid) stub per layer per track.
0786         const StubStreamData& stubdata = streamsStubRaw[chanStubOffsetIn + iproj][itk];
0787         const L1TStub& stub = stubdata.stub();
0788         if (!stubdata.valid())
0789           continue;
0790         const TTStubRef& ttStubRef = stubMap[stub];
0791         const int seedType = stubdata.iSeed();
0792         const int layerId = setup_->layerId(ttStubRef);
0793         const int channelId = channelAssignment_->channelId(seedType, layerId);
0794         hitMap.set(channelId);
0795         streamsStub[chanStubOffsetOut + channelId].emplace_back(ttStubRef, stubdata.dataBits());
0796       }
0797       for (int layerId : hitMap.ids(false)) {  // invalid stubs
0798         streamsStub[chanStubOffsetOut + layerId].emplace_back(tt::FrameStub());
0799       }
0800     }
0801   }
0802 
0803   iEvent.emplace(putTokenTracks_, std::move(streamsTrack));
0804   iEvent.emplace(putTokenStubs_, std::move(streamsStub));
0805 
0806 }  /// End of produce()
0807 
0808 // ///////////////////////////
0809 // // DEFINE THIS AS A PLUG-IN
0810 DEFINE_FWK_MODULE(L1FPGATrackProducer);