Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #ifndef L1Trigger_TrackFindingTMTT_Stub_h
0002 #define L1Trigger_TrackFindingTMTT_Stub_h
0003 
0004 #include "DataFormats/L1TrackTrigger/interface/TTTypes.h"
0005 #include "DataFormats/L1TrackTrigger/interface/TTStub.h"
0006 #include "DataFormats/Phase2TrackerDigi/interface/Phase2TrackerDigi.h"
0007 #include "SimDataFormats/Associations/interface/TTStubAssociationMap.h"
0008 #include "DataFormats/DetId/interface/DetId.h"
0009 #include "DataFormats/Common/interface/Ref.h"
0010 #include "DataFormats/Common/interface/DetSetVector.h"
0011 #include "FWCore/Framework/interface/Frameworkfwd.h"
0012 
0013 #include "L1Trigger/TrackFindingTMTT/interface/Settings.h"
0014 #include "L1Trigger/TrackFindingTMTT/interface/DigitalStub.h"
0015 #include "L1Trigger/TrackFindingTMTT/interface/DegradeBend.h"
0016 #include "L1Trigger/TrackFindingTMTT/interface/TrackerModule.h"
0017 
0018 #include <vector>
0019 #include <set>
0020 #include <array>
0021 #include <map>
0022 #include <memory>
0023 
0024 class TrackerGeometry;
0025 class TrackerTopology;
0026 
0027 namespace tmtt {
0028 
0029   class TP;
0030   class DegradeBend;
0031   class StubKiller;
0032 
0033   typedef edmNew::DetSetVector<TTStub<Ref_Phase2TrackerDigi_> > TTStubDetSetVec;
0034   typedef edmNew::DetSet<TTStub<Ref_Phase2TrackerDigi_> > TTStubDetSet;
0035   typedef edm::Ref<TTStubDetSetVec, TTStub<Ref_Phase2TrackerDigi_> > TTStubRef;
0036   typedef edm::Ref<edmNew::DetSetVector<TTCluster<Ref_Phase2TrackerDigi_> >, TTCluster<Ref_Phase2TrackerDigi_> >
0037       TTClusterRef;
0038   typedef TTStubAssociationMap<Ref_Phase2TrackerDigi_> TTStubAssMap;
0039   typedef TTClusterAssociationMap<Ref_Phase2TrackerDigi_> TTClusterAssMap;
0040 
0041   //=== Represents a Tracker stub (=pair of hits)
0042 
0043   class Stub {
0044   public:
0045     // Hybrid L1 tracking: stub constructor.
0046     Stub(const Settings* settings,
0047          unsigned int idStub,
0048          double phi,
0049          double r,
0050          double z,
0051          double bend,
0052          unsigned int iphi,
0053          double alpha,
0054          unsigned int layerId,
0055          unsigned int iPhiSec,
0056          bool psModule,
0057          bool barrel,
0058          bool tiltedBarrel,
0059          float stripPitch,
0060          float stripLength,
0061          unsigned int nStrips);
0062 
0063     // TMTT L1 tracking: stub constructor.
0064     Stub(const TTStubRef& ttStubRef,
0065          unsigned int index_in_vStubs,
0066          const Settings* settings,
0067          const TrackerTopology* trackerTopology,
0068          const TrackerModule* trackerModule,
0069          const DegradeBend* degradeBend,
0070          const StubKiller* stubKiller);
0071 
0072     // KF emualtor: stub constructor
0073     Stub(const TTStubRef& ttStubRef,
0074          double r,
0075          double phi,
0076          double z,
0077          int layerId,
0078          int layerIdReduced,
0079          double stripPitch,
0080          double stripLength,
0081          bool psModule,
0082          bool barrel,
0083          bool tiltedBarrel);
0084 
0085     bool operator==(const Stub& stubOther) { return (this->index() == stubOther.index()); }
0086 
0087     // Return reference to original TTStub.
0088     const TTStubRef& ttStubRef() const { return ttStubRef_; }
0089 
0090     // Info about tracker module containing stub.
0091     const TrackerModule* trackerModule() const { return trackerModule_; }
0092 
0093     // Fill truth info with association from stub to tracking particles.
0094     void fillTruth(const std::map<edm::Ptr<TrackingParticle>, const TP*>& translateTP,
0095                    const edm::Handle<TTStubAssMap>& mcTruthTTStubHandle,
0096                    const edm::Handle<TTClusterAssMap>& mcTruthTTClusterHandle);
0097 
0098     // Calculate HT m-bin range consistent with bend.
0099     void calcQoverPtrange();
0100 
0101     // Digitize stub for input to GP, HT, SF, TF
0102     enum class DigiStage { NONE, GP, HT, SF, TF };
0103     void digitize(unsigned int iPhiSec, DigiStage digiStep);
0104 
0105     // Control warning messages about accessing non-digitized quantities.
0106     void setDigitizeWarningsOn(bool newVal) { digitizeWarningsOn_ = newVal; }
0107 
0108     // Access to digitized version of stub coords.
0109     const DigitalStub* digitalStub() const { return digitalStub_.get(); }
0110 
0111     // === Functions for returning info about reconstructed stubs ===
0112 
0113     // Location in InputData::vStubs_
0114     unsigned int index() const { return index_in_vStubs_; }
0115 
0116     //--- Stub data and quantities derived from it ---
0117 
0118     // Stub coordinates (optionally after digitisation, if digitisation requested via cfg).
0119     // N.B. Digitisation is only run if Stub::digitize() is called.
0120     float phi() const { return phi_; }
0121     float r() const { return r_; }
0122     float z() const { return z_; }
0123     float theta() const { return atan2(r_, z_); }
0124     float eta() const { return asinh(z_ / r_); }
0125 
0126     // Location of stub in module in units of strip/pixel number in phi direction.
0127     // Range from 0 to (nStrips - 1) inclusive.
0128     unsigned int iphi() const { return iphi_; }
0129     // alpha correction for non-radial strips in endcap 2S modules.
0130     // (If true hit at larger r than stub r by deltaR, then stub phi needs correcting by +alpha*deltaR).
0131     // *** TO DO *** : Digitize this.
0132     float alpha() const { return alpha_; }
0133 
0134     // Get stub bend and its resolution, as available within the front end chip (i.e. prior to loss of bits
0135     // or digitisation).
0136     float bendInFrontend() const { return bendInFrontend_; }
0137     float bendCutInFrontend() const { return settings_->bendCut(); }
0138     // Get stub bend (i.e. displacement between two hits in stub in units of strip pitch).
0139     float bend() const { return bend_; }
0140     // Bend resolution.
0141     float bendCut() const { return (settings_->bendCut() + (numMergedBend_ - 1) * settings_->bendCutExtra()); }
0142     // No. of bend values merged into FE bend encoding of this stub.
0143     float numMergedBend() const { return numMergedBend_; }
0144     // Estimated track q/Pt based on stub bend info.
0145     float qOverPt() const { return (this->qOverPtOverBend() * this->bend()); }
0146     float qOverPtcut() const { return (this->qOverPtOverBend() * this->bendCut()); }
0147     // Range in q/Pt bins in HT array compatible with stub bend.
0148     unsigned int min_qOverPt_bin() const { return min_qOverPt_bin_; }
0149     unsigned int max_qOverPt_bin() const { return max_qOverPt_bin_; }
0150     // Difference in phi between stub and angle at which track crosses given radius, assuming track has given Pt.
0151     float phiDiff(float rad, float Pt) const { return std::abs(r_ - rad) * (settings_->invPtToDphi()) / Pt; }
0152     // Phi angle at which particle consistent with this stub & its bend cross specified radius.
0153     float trkPhiAtR(float rad) const { return phi_ + (bend_ * dphiOverBend_) * (1. - rad / r_); }
0154     // Its resolution
0155     float trkPhiAtRcut(float rad) const { return (bendCut() * dphiOverBend_) * std::abs(1. - rad / r_); }
0156 
0157     // -- conversion factors
0158     // Ratio of track crossing angle to bend.
0159     float dphiOverBend() const { return dphiOverBend_; }
0160     // Ratio of q/Pt to bend.
0161     float qOverPtOverBend() const { return dphiOverBend_ / (r_ * settings_->invPtToDphi()); }
0162 
0163     //--- Info about the two clusters that make up the stub.
0164 
0165     // Coordinates in frame of sensor, measured in units of strip pitch along two orthogonal axes running perpendicular and parallel to longer axis of pixels/strips (U & V).
0166     std::array<float, 2> localU_cluster() const { return localU_cluster_; }
0167     std::array<float, 2> localV_cluster() const { return localV_cluster_; }
0168 
0169     //--- Check if this stub will be output by FE. Stub failing this not used for L1 tracks.
0170     bool frontendPass() const { return frontendPass_; }
0171     // Indicates if stub would have passed DE cuts, were it not for window size encoded in DegradeBend.h
0172     bool stubFailedDegradeWindow() const { return stubFailedDegradeWindow_; }
0173 
0174     //--- Truth info
0175 
0176     // Association of stub to tracking particles
0177     const std::set<const TP*>& assocTPs() const {
0178       return assocTPs_;
0179     }  // Return TPs associated to this stub. (Whether only TPs contributing to both clusters are returned is determined by "StubMatchStrict" config param.)
0180     bool genuine() const { return (not assocTPs_.empty()); }  // Did stub match at least one TP?
0181     const TP* assocTP() const {
0182       return assocTP_;
0183     }  // If only one TP contributed to both clusters, this tells you which TP it is. Returns nullptr if none.
0184 
0185     // Association of both clusters making up stub to tracking particles
0186     std::array<bool, 2> genuineCluster() const {
0187       return std::array<bool, 2>{{(assocTPofCluster_[0] != nullptr), (assocTPofCluster_[1] != nullptr)}};
0188     }  // Was cluster produced by a single TP?
0189     std::array<const TP*, 2> assocTPofCluster() const {
0190       return assocTPofCluster_;
0191     }  // Which TP made each cluster. Warning: If cluster was not produced by a single TP, then returns nullptr! (P.S. If both clusters match same TP, then this will equal assocTP()).
0192 
0193     //--- Quantities common to all stubs in a given module ---
0194     // N.B. Not taken from trackerModule_ to cope with Hybrid tracking.
0195 
0196     // Angle between normal to module and beam-line along +ve z axis. (In range -PI/2 to +PI/2).
0197     float tiltAngle() const { return tiltAngle_; }
0198     // Uncertainty in stub (r,z)
0199     float sigmaR() const { return (barrel() ? 0. : sigmaPar()); }
0200     float sigmaZ() const { return (barrel() ? sigmaPar() : 0.); }
0201     // Hit resolution perpendicular to strip. Measures phi.
0202     float sigmaPerp() const { return invRoot12 * stripPitch_; }
0203     // Hit resolution parallel to strip. Measures r or z.
0204     float sigmaPar() const { return invRoot12 * stripLength_; }
0205 
0206     //--- These module variables could be taken directly from trackerModule_, were it not for need
0207     //--- to support Hybrid.
0208     // Module type: PS or 2S?
0209     bool psModule() const { return psModule_; }
0210     // Tracker layer ID number (1-6 = barrel layer; 11-15 = endcap A disk; 21-25 = endcap B disk)
0211     unsigned int layerId() const { return layerId_; }
0212     // Reduced layer ID (in range 1-7). This encodes the layer ID in only 3 bits (to simplify firmware) by merging some barrel layer and endcap disk layer IDs into a single ID.
0213     unsigned int layerIdReduced() const { return layerIdReduced_; }
0214     bool barrel() const { return barrel_; }
0215     // True if stub is in tilted barrel module.
0216     bool tiltedBarrel() const { return tiltedBarrel_; }
0217     // Strip pitch (or pixel pitch along shortest axis).
0218     float stripPitch() const { return stripPitch_; }
0219     // Strip length (or pixel pitch along longest axis).
0220     float stripLength() const { return stripLength_; }
0221     // No. of strips in sensor.
0222     unsigned int nStrips() const { return nStrips_; }
0223 
0224   private:
0225     // Degrade assumed stub bend resolution.
0226     // And return an integer indicating how many values of bend are merged into this single one.
0227     void degradeResolution(float bend, float& degradedBend, unsigned int& num) const;
0228 
0229     // Set the frontendPass_ flag, indicating if frontend readout electronics will output this stub.
0230     void setFrontend(const StubKiller* stubKiller);
0231 
0232     // Set info about the module that this stub is in.
0233     void setTrackerModule(const TrackerGeometry* trackerGeometry,
0234                           const TrackerTopology* trackerTopology,
0235                           const DetId& detId);
0236 
0237     // Function to calculate approximation for dphiOverBendCorrection aka B
0238     double approxB();
0239 
0240     // Calculate variables giving ratio of track intercept angle to stub bend.
0241     void calcDphiOverBend();
0242 
0243   private:
0244     TTStubRef ttStubRef_;  // Reference to original TTStub
0245 
0246     const Settings* settings_;  // configuration parameters.
0247 
0248     unsigned int index_in_vStubs_;  // location of this stub in InputData::vStubs
0249 
0250     //--- Parameters passed along optical links from PP to MP (or equivalent ones if easier for analysis software to use).
0251     // WARNING: If you add any variables in this section, take care to ensure that they are digitized correctly by Stub::digitize().
0252     float phi_;  // stub coords, optionally after digitisation.
0253     float r_;
0254     float z_;
0255     float bend_;                    // bend of stub.
0256     float dphiOverBend_;            // related to rho parameter.
0257     unsigned int min_qOverPt_bin_;  // Range in q/Pt bins in HT array compatible with stub bend.
0258     unsigned int max_qOverPt_bin_;
0259 
0260     //--- Info about the two clusters that make up the stub.
0261     std::array<float, 2> localU_cluster_;
0262     std::array<float, 2> localV_cluster_;
0263 
0264     unsigned int iphi_;
0265     float alpha_;
0266 
0267     // Would front-end electronics output this stub?
0268     bool frontendPass_;
0269     // Did stub fail window cuts assumed in DegradeBend.h?
0270     bool stubFailedDegradeWindow_;
0271     // Bend in front end chip (prior to degredation by loss of bits & digitization).
0272     float bendInFrontend_;
0273     // Used for stub bend resolution degrading.
0274     unsigned int numMergedBend_;
0275 
0276     //--- Truth info about stub.
0277     const TP* assocTP_;
0278     std::set<const TP*> assocTPs_;
0279     //--- Truth info about the two clusters that make up the stub
0280     std::array<const TP*, 2> assocTPofCluster_;
0281 
0282     std::unique_ptr<DigitalStub> digitalStub_;  // Class used to digitize stub if required.
0283     DigiStage lastDigiStep_;
0284     bool digitizeWarningsOn_;  // Enable warnings about accessing non-digitized quantities.
0285 
0286     // Info about tracker module containing stub.
0287     const TrackerModule* trackerModule_;
0288 
0289     // Used to degrade stub bend information.
0290     const DegradeBend* degradeBend_;
0291 
0292     // These module variables are needed only to support the Hybrid stub constructor.
0293     // (Otherwise, they could be taken from trackerModule_).
0294     unsigned int layerId_;
0295     unsigned int layerIdReduced_;
0296     float tiltAngle_;
0297     float stripPitch_;
0298     float stripLength_;
0299     unsigned int nStrips_;
0300     bool psModule_;
0301     bool barrel_;
0302     bool tiltedBarrel_;
0303 
0304     const float rejectedStubBend_ = 99999.;  // Bend set to this if stub rejected.
0305 
0306     const float invRoot12 = sqrt(1. / 12.);
0307   };
0308 
0309 }  // namespace tmtt
0310 #endif