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
0042
0043 class Stub {
0044 public:
0045
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
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
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
0088 const TTStubRef& ttStubRef() const { return ttStubRef_; }
0089
0090
0091 const TrackerModule* trackerModule() const { return trackerModule_; }
0092
0093
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
0099 void calcQoverPtrange();
0100
0101
0102 enum class DigiStage { NONE, GP, HT, SF, TF };
0103 void digitize(unsigned int iPhiSec, DigiStage digiStep);
0104
0105
0106 void setDigitizeWarningsOn(bool newVal) { digitizeWarningsOn_ = newVal; }
0107
0108
0109 const DigitalStub* digitalStub() const { return digitalStub_.get(); }
0110
0111
0112
0113
0114 unsigned int index() const { return index_in_vStubs_; }
0115
0116
0117
0118
0119
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
0127
0128 unsigned int iphi() const { return iphi_; }
0129
0130
0131
0132 float alpha() const { return alpha_; }
0133
0134
0135
0136 float bendInFrontend() const { return bendInFrontend_; }
0137 float bendCutInFrontend() const { return settings_->bendCut(); }
0138
0139 float bend() const { return bend_; }
0140
0141 float bendCut() const { return (settings_->bendCut() + (numMergedBend_ - 1) * settings_->bendCutExtra()); }
0142
0143 float numMergedBend() const { return numMergedBend_; }
0144
0145 float qOverPt() const { return (this->qOverPtOverBend() * this->bend()); }
0146 float qOverPtcut() const { return (this->qOverPtOverBend() * this->bendCut()); }
0147
0148 unsigned int min_qOverPt_bin() const { return min_qOverPt_bin_; }
0149 unsigned int max_qOverPt_bin() const { return max_qOverPt_bin_; }
0150
0151 float phiDiff(float rad, float Pt) const { return std::abs(r_ - rad) * (settings_->invPtToDphi()) / Pt; }
0152
0153 float trkPhiAtR(float rad) const { return phi_ + (bend_ * dphiOverBend_) * (1. - rad / r_); }
0154
0155 float trkPhiAtRcut(float rad) const { return (bendCut() * dphiOverBend_) * std::abs(1. - rad / r_); }
0156
0157
0158
0159 float dphiOverBend() const { return dphiOverBend_; }
0160
0161 float qOverPtOverBend() const { return dphiOverBend_ / (r_ * settings_->invPtToDphi()); }
0162
0163
0164
0165
0166 std::array<float, 2> localU_cluster() const { return localU_cluster_; }
0167 std::array<float, 2> localV_cluster() const { return localV_cluster_; }
0168
0169
0170 bool frontendPass() const { return frontendPass_; }
0171
0172 bool stubFailedDegradeWindow() const { return stubFailedDegradeWindow_; }
0173
0174
0175
0176
0177 const std::set<const TP*>& assocTPs() const {
0178 return assocTPs_;
0179 }
0180 bool genuine() const { return (not assocTPs_.empty()); }
0181 const TP* assocTP() const {
0182 return assocTP_;
0183 }
0184
0185
0186 std::array<bool, 2> genuineCluster() const {
0187 return std::array<bool, 2>{{(assocTPofCluster_[0] != nullptr), (assocTPofCluster_[1] != nullptr)}};
0188 }
0189 std::array<const TP*, 2> assocTPofCluster() const {
0190 return assocTPofCluster_;
0191 }
0192
0193
0194
0195
0196
0197 float tiltAngle() const { return tiltAngle_; }
0198
0199 float sigmaR() const { return (barrel() ? 0. : sigmaPar()); }
0200 float sigmaZ() const { return (barrel() ? sigmaPar() : 0.); }
0201
0202 float sigmaPerp() const { return invRoot12 * stripPitch_; }
0203
0204 float sigmaPar() const { return invRoot12 * stripLength_; }
0205
0206
0207
0208
0209 bool psModule() const { return psModule_; }
0210
0211 unsigned int layerId() const { return layerId_; }
0212
0213 unsigned int layerIdReduced() const { return layerIdReduced_; }
0214 bool barrel() const { return barrel_; }
0215
0216 bool tiltedBarrel() const { return tiltedBarrel_; }
0217
0218 float stripPitch() const { return stripPitch_; }
0219
0220 float stripLength() const { return stripLength_; }
0221
0222 unsigned int nStrips() const { return nStrips_; }
0223
0224 private:
0225
0226
0227 void degradeResolution(float bend, float& degradedBend, unsigned int& num) const;
0228
0229
0230 void setFrontend(const StubKiller* stubKiller);
0231
0232
0233 void setTrackerModule(const TrackerGeometry* trackerGeometry,
0234 const TrackerTopology* trackerTopology,
0235 const DetId& detId);
0236
0237
0238 double approxB();
0239
0240
0241 void calcDphiOverBend();
0242
0243 private:
0244 TTStubRef ttStubRef_;
0245
0246 const Settings* settings_;
0247
0248 unsigned int index_in_vStubs_;
0249
0250
0251
0252 float phi_;
0253 float r_;
0254 float z_;
0255 float bend_;
0256 float dphiOverBend_;
0257 unsigned int min_qOverPt_bin_;
0258 unsigned int max_qOverPt_bin_;
0259
0260
0261 std::array<float, 2> localU_cluster_;
0262 std::array<float, 2> localV_cluster_;
0263
0264 unsigned int iphi_;
0265 float alpha_;
0266
0267
0268 bool frontendPass_;
0269
0270 bool stubFailedDegradeWindow_;
0271
0272 float bendInFrontend_;
0273
0274 unsigned int numMergedBend_;
0275
0276
0277 const TP* assocTP_;
0278 std::set<const TP*> assocTPs_;
0279
0280 std::array<const TP*, 2> assocTPofCluster_;
0281
0282 std::unique_ptr<DigitalStub> digitalStub_;
0283 DigiStage lastDigiStep_;
0284 bool digitizeWarningsOn_;
0285
0286
0287 const TrackerModule* trackerModule_;
0288
0289
0290 const DegradeBend* degradeBend_;
0291
0292
0293
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.;
0305
0306 const float invRoot12 = sqrt(1. / 12.);
0307 };
0308
0309 }
0310 #endif