Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-06-03 00:12:21

0001 #ifndef L1Trigger_TrackFindingTracklet_TrackMultiplexer_h
0002 #define L1Trigger_TrackFindingTracklet_TrackMultiplexer_h
0003 
0004 #include "L1Trigger/TrackTrigger/interface/Setup.h"
0005 #include "L1Trigger/TrackFindingTracklet/interface/ChannelAssignment.h"
0006 #include "L1Trigger/TrackFindingTracklet/interface/DataFormats.h"
0007 #include "L1Trigger/TrackFindingTracklet/interface/Settings.h"
0008 
0009 #include <vector>
0010 #include <deque>
0011 
0012 namespace trklet {
0013 
0014   /*! \class  trklet::TrackMultiplexer
0015    *  \brief  Class to emulate transformation of tracklet tracks and stubs into TMTT format
0016    *          and routing of seed type streams into single stream
0017    *  \author Thomas Schuh
0018    *  \date   2023, Jan
0019    */
0020   class TrackMultiplexer {
0021   public:
0022     TrackMultiplexer(const tt::Setup* setup_,
0023                      const DataFormats* dataFormats,
0024                      const ChannelAssignment* channelAssignment,
0025                      const Settings* settings,
0026                      int region);
0027     ~TrackMultiplexer() = default;
0028     // read in and organize input tracks and stubs
0029     void consume(const tt::StreamsTrack& streamsTrack, const tt::StreamsStub& streamsStub);
0030     // fill output products
0031     void produce(tt::StreamsTrack& streamsTrack, tt::StreamsStub& streamsStub);
0032 
0033   private:
0034     // truncates double precision of val into base precision, +1.e-12 restores robustness of addition of 2 digitised values
0035     double digi(double val, double base) const { return (floor(val / base + 1.e-12) + .5) * base; }
0036     // basetransformation of val from baseLow into baseHigh using widthMultiplier bit multiplication
0037     double redigi(double val, double baseLow, double baseHigh, int widthMultiplier) const;
0038     struct Stub {
0039       Stub(const TTStubRef& ttStubRef, int layer, int stubId, double r, double phi, double z, bool psTilt)
0040           : valid_(true), ttStubRef_(ttStubRef), layer_(layer), stubId_(stubId), r_(r), phi_(phi), z_(z) {
0041         stubId_ = 2 * stubId_ + (psTilt ? 1 : 0);
0042       }
0043       tt::FrameStub frame(const DataFormats* df) const { return StubTM(ttStubRef_, df, stubId_, r_, phi_, z_).frame(); }
0044       bool valid_;
0045       TTStubRef ttStubRef_;
0046       // kf layer id
0047       int layer_;
0048       // tracklet stub id, used to identify duplicates
0049       int stubId_;
0050       // radius w.r.t. chosenRofPhi in cm
0051       double r_;
0052       // phi residual in rad
0053       double phi_;
0054       // z residual in cm
0055       double z_;
0056     };
0057     struct Track {
0058       static constexpr int max_ = 11;
0059       Track() { stubs_.reserve(max_); }
0060       Track(const TTTrackRef& ttTrackRef,
0061             bool valid,
0062             int seedType,
0063             double inv2R,
0064             double phiT,
0065             double cot,
0066             double zT,
0067             const std::vector<Stub*>& stubs)
0068           : ttTrackRef_(ttTrackRef),
0069             valid_(valid),
0070             seedType_(seedType),
0071             inv2R_(inv2R),
0072             phiT_(phiT),
0073             cot_(cot),
0074             zT_(zT),
0075             stubs_(stubs) {}
0076       tt::FrameTrack frame(const DataFormats* df) const { return TrackTM(ttTrackRef_, df, inv2R_, phiT_, zT_).frame(); }
0077       TTTrackRef ttTrackRef_;
0078       bool valid_;
0079       int seedType_;
0080       double inv2R_;
0081       double phiT_;
0082       double cot_;
0083       double zT_;
0084       std::vector<Stub*> stubs_;
0085     };
0086     // remove and return first element of deque, returns nullptr if empty
0087     template <class T>
0088     T* pop_front(std::deque<T*>& ts) const;
0089     // true if truncation is enbaled
0090     bool enableTruncation_;
0091     // stub residuals are recalculated from seed parameter and TTStub position
0092     bool useTTStubResiduals_;
0093     // track parameter are recalculated from seed TTStub positions
0094     bool useTTStubParameters_;
0095     //
0096     bool applyNonLinearCorrection_;
0097     // provides run-time constants
0098     const tt::Setup* setup_;
0099     // provides dataformats
0100     const DataFormats* dataFormats_;
0101     // helper class to assign tracks to channel
0102     const ChannelAssignment* channelAssignment_;
0103     // provides tracklet constants
0104     const Settings* settings_;
0105     // processing region (0 - 8) aka processing phi nonant
0106     const int region_;
0107     // storage of input tracks
0108     std::vector<Track> tracks_;
0109     // storage of input stubs
0110     std::vector<Stub> stubs_;
0111     // h/w liked organized pointer to input tracks
0112     std::vector<std::vector<Track*>> input_;
0113     // unified tracklet digitisation granularity
0114     double baseUinv2R_;
0115     double baseUphiT_;
0116     double baseUcot_;
0117     double baseUzT_;
0118     double baseUr_;
0119     double baseUphi_;
0120     double baseUz_;
0121     // KF input format digitisation granularity (identical to TMTT)
0122     double baseLinv2R_;
0123     double baseLphiT_;
0124     double baseLzT_;
0125     double baseLr_;
0126     double baseLphi_;
0127     double baseLz_;
0128     double baseLcot_;
0129     // Finer granularity (by powers of 2) than the TMTT one. Used to transform from Tracklet to TMTT base.
0130     double baseHinv2R_;
0131     double baseHphiT_;
0132     double baseHzT_;
0133     double baseHr_;
0134     double baseHphi_;
0135     double baseHz_;
0136     double baseHcot_;
0137     // digitisation granularity used for inverted cot(theta)
0138     double baseInvCot_;
0139     double baseScot_;
0140   };
0141 
0142 }  // namespace trklet
0143 
0144 #endif