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