CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Stub.h
Go to the documentation of this file.
1 #ifndef L1Trigger_TrackFindingTMTT_Stub_h
2 #define L1Trigger_TrackFindingTMTT_Stub_h
3 
12 
17 
18 #include <vector>
19 #include <set>
20 #include <array>
21 #include <map>
22 #include <memory>
23 
24 class TrackerGeometry;
25 class TrackerTopology;
26 
27 namespace tmtt {
28 
29  class TP;
30  class DegradeBend;
31  class StubKiller;
32 
40 
41  //=== Represents a Tracker stub (=pair of hits)
42 
43  class Stub {
44  public:
45  // Hybrid L1 tracking: stub constructor.
46  Stub(const Settings* settings,
47  unsigned int idStub,
48  double phi,
49  double r,
50  double z,
51  double bend,
52  unsigned int iphi,
53  double alpha,
54  unsigned int layerId,
55  unsigned int iPhiSec,
56  bool psModule,
57  bool barrel,
58  bool tiltedBarrel,
59  float stripPitch,
60  float stripLength,
61  unsigned int nStrips);
62 
63  // TMTT L1 tracking: stub constructor.
64  Stub(const TTStubRef& ttStubRef,
65  unsigned int index_in_vStubs,
66  const Settings* settings,
67  const TrackerTopology* trackerTopology,
69  const DegradeBend* degradeBend,
70  const StubKiller* stubKiller);
71 
72  bool operator==(const Stub& stubOther) { return (this->index() == stubOther.index()); }
73 
74  // Return reference to original TTStub.
75  const TTStubRef& ttStubRef() const { return ttStubRef_; }
76 
77  // Info about tracker module containing stub.
78  const TrackerModule* trackerModule() const { return trackerModule_; }
79 
80  // Fill truth info with association from stub to tracking particles.
81  void fillTruth(const std::map<edm::Ptr<TrackingParticle>, const TP*>& translateTP,
82  const edm::Handle<TTStubAssMap>& mcTruthTTStubHandle,
83  const edm::Handle<TTClusterAssMap>& mcTruthTTClusterHandle);
84 
85  // Calculate HT m-bin range consistent with bend.
86  void calcQoverPtrange();
87 
88  // Digitize stub for input to GP, HT, SF, TF
89  enum class DigiStage { NONE, GP, HT, SF, TF };
90  void digitize(unsigned int iPhiSec, DigiStage digiStep);
91 
92  // Control warning messages about accessing non-digitized quantities.
93  void setDigitizeWarningsOn(bool newVal) { digitizeWarningsOn_ = newVal; }
94 
95  // Access to digitized version of stub coords.
96  const DigitalStub* digitalStub() const { return digitalStub_.get(); }
97 
98  // === Functions for returning info about reconstructed stubs ===
99 
100  // Location in InputData::vStubs_
101  unsigned int index() const { return index_in_vStubs_; }
102 
103  //--- Stub data and quantities derived from it ---
104 
105  // Stub coordinates (optionally after digitisation, if digitisation requested via cfg).
106  // N.B. Digitisation is only run if Stub::digitize() is called.
107  float phi() const { return phi_; }
108  float r() const { return r_; }
109  float z() const { return z_; }
110  float theta() const { return atan2(r_, z_); }
111  float eta() const { return asinh(z_ / r_); }
112 
113  // Location of stub in module in units of strip/pixel number in phi direction.
114  // Range from 0 to (nStrips - 1) inclusive.
115  unsigned int iphi() const { return iphi_; }
116  // alpha correction for non-radial strips in endcap 2S modules.
117  // (If true hit at larger r than stub r by deltaR, then stub phi needs correcting by +alpha*deltaR).
118  // *** TO DO *** : Digitize this.
119  float alpha() const { return alpha_; }
120 
121  // Get stub bend and its resolution, as available within the front end chip (i.e. prior to loss of bits
122  // or digitisation).
123  float bendInFrontend() const { return bendInFrontend_; }
124  float bendCutInFrontend() const { return settings_->bendCut(); }
125  // Get stub bend (i.e. displacement between two hits in stub in units of strip pitch).
126  float bend() const { return bend_; }
127  // Bend resolution.
128  float bendCut() const { return (settings_->bendCut() + (numMergedBend_ - 1) * settings_->bendCutExtra()); }
129  // No. of bend values merged into FE bend encoding of this stub.
130  float numMergedBend() const { return numMergedBend_; }
131  // Estimated track q/Pt based on stub bend info.
132  float qOverPt() const { return (this->qOverPtOverBend() * this->bend()); }
133  float qOverPtcut() const { return (this->qOverPtOverBend() * this->bendCut()); }
134  // Range in q/Pt bins in HT array compatible with stub bend.
135  unsigned int min_qOverPt_bin() const { return min_qOverPt_bin_; }
136  unsigned int max_qOverPt_bin() const { return max_qOverPt_bin_; }
137  // Difference in phi between stub and angle at which track crosses given radius, assuming track has given Pt.
138  float phiDiff(float rad, float Pt) const { return std::abs(r_ - rad) * (settings_->invPtToDphi()) / Pt; }
139  // Phi angle at which particle consistent with this stub & its bend cross specified radius.
140  float trkPhiAtR(float rad) const { return phi_ + (bend_ * dphiOverBend_) * (1. - rad / r_); }
141  // Its resolution
142  float trkPhiAtRcut(float rad) const { return (bendCut() * dphiOverBend_) * std::abs(1. - rad / r_); }
143 
144  // -- conversion factors
145  // Ratio of track crossing angle to bend.
146  float dphiOverBend() const { return dphiOverBend_; }
147  // Ratio of q/Pt to bend.
148  float qOverPtOverBend() const { return dphiOverBend_ / (r_ * settings_->invPtToDphi()); }
149 
150  //--- Info about the two clusters that make up the stub.
151 
152  // 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).
153  std::array<float, 2> localU_cluster() const { return localU_cluster_; }
154  std::array<float, 2> localV_cluster() const { return localV_cluster_; }
155 
156  //--- Check if this stub will be output by FE. Stub failing this not used for L1 tracks.
157  bool frontendPass() const { return frontendPass_; }
158  // Indicates if stub would have passed DE cuts, were it not for window size encoded in DegradeBend.h
160 
161  //--- Truth info
162 
163  // Association of stub to tracking particles
164  const std::set<const TP*>& assocTPs() const {
165  return assocTPs_;
166  } // Return TPs associated to this stub. (Whether only TPs contributing to both clusters are returned is determined by "StubMatchStrict" config param.)
167  bool genuine() const { return (not assocTPs_.empty()); } // Did stub match at least one TP?
168  const TP* assocTP() const {
169  return assocTP_;
170  } // If only one TP contributed to both clusters, this tells you which TP it is. Returns nullptr if none.
171 
172  // Association of both clusters making up stub to tracking particles
173  std::array<bool, 2> genuineCluster() const {
174  return std::array<bool, 2>{{(assocTPofCluster_[0] != nullptr), (assocTPofCluster_[1] != nullptr)}};
175  } // Was cluster produced by a single TP?
176  std::array<const TP*, 2> assocTPofCluster() const {
177  return assocTPofCluster_;
178  } // 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()).
179 
180  //--- Quantities common to all stubs in a given module ---
181  // N.B. Not taken from trackerModule_ to cope with Hybrid tracking.
182 
183  // Angle between normal to module and beam-line along +ve z axis. (In range -PI/2 to +PI/2).
184  float tiltAngle() const { return tiltAngle_; }
185  // Uncertainty in stub (r,z)
186  float sigmaR() const { return (barrel() ? 0. : sigmaPar()); }
187  float sigmaZ() const { return (barrel() ? sigmaPar() : 0.); }
188  // Hit resolution perpendicular to strip. Measures phi.
189  float sigmaPerp() const { return invRoot12 * stripPitch_; }
190  // Hit resolution parallel to strip. Measures r or z.
191  float sigmaPar() const { return invRoot12 * stripLength_; }
192 
193  //--- These module variables could be taken directly from trackerModule_, were it not for need
194  //--- to support Hybrid.
195  // Module type: PS or 2S?
196  bool psModule() const { return psModule_; }
197  // Tracker layer ID number (1-6 = barrel layer; 11-15 = endcap A disk; 21-25 = endcap B disk)
198  unsigned int layerId() const { return layerId_; }
199  // 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.
200  unsigned int layerIdReduced() const { return layerIdReduced_; }
201  bool barrel() const { return barrel_; }
202  // True if stub is in tilted barrel module.
203  bool tiltedBarrel() const { return tiltedBarrel_; }
204  // Strip pitch (or pixel pitch along shortest axis).
205  float stripPitch() const { return stripPitch_; }
206  // Strip length (or pixel pitch along longest axis).
207  float stripLength() const { return stripLength_; }
208  // No. of strips in sensor.
209  unsigned int nStrips() const { return nStrips_; }
210 
211  private:
212  // Degrade assumed stub bend resolution.
213  // And return an integer indicating how many values of bend are merged into this single one.
214  void degradeResolution(float bend, float& degradedBend, unsigned int& num) const;
215 
216  // Set the frontendPass_ flag, indicating if frontend readout electronics will output this stub.
217  void setFrontend(const StubKiller* stubKiller);
218 
219  // Set info about the module that this stub is in.
220  void setTrackerModule(const TrackerGeometry* trackerGeometry,
221  const TrackerTopology* trackerTopology,
222  const DetId& detId);
223 
224  // Function to calculate approximation for dphiOverBendCorrection aka B
225  double approxB();
226 
227  // Calculate variables giving ratio of track intercept angle to stub bend.
228  void calcDphiOverBend();
229 
230  private:
231  TTStubRef ttStubRef_; // Reference to original TTStub
232 
233  const Settings* settings_; // configuration parameters.
234 
235  unsigned int index_in_vStubs_; // location of this stub in InputData::vStubs
236 
237  //--- Parameters passed along optical links from PP to MP (or equivalent ones if easier for analysis software to use).
238  // WARNING: If you add any variables in this section, take care to ensure that they are digitized correctly by Stub::digitize().
239  float phi_; // stub coords, optionally after digitisation.
240  float r_;
241  float z_;
242  float bend_; // bend of stub.
243  float dphiOverBend_; // related to rho parameter.
244  unsigned int min_qOverPt_bin_; // Range in q/Pt bins in HT array compatible with stub bend.
245  unsigned int max_qOverPt_bin_;
246 
247  //--- Info about the two clusters that make up the stub.
248  std::array<float, 2> localU_cluster_;
249  std::array<float, 2> localV_cluster_;
250 
251  unsigned int iphi_;
252  float alpha_;
253 
254  // Would front-end electronics output this stub?
256  // Did stub fail window cuts assumed in DegradeBend.h?
258  // Bend in front end chip (prior to degredation by loss of bits & digitization).
260  // Used for stub bend resolution degrading.
261  unsigned int numMergedBend_;
262 
263  //--- Truth info about stub.
264  const TP* assocTP_;
265  std::set<const TP*> assocTPs_;
266  //--- Truth info about the two clusters that make up the stub
267  std::array<const TP*, 2> assocTPofCluster_;
268 
269  std::unique_ptr<DigitalStub> digitalStub_; // Class used to digitize stub if required.
271  bool digitizeWarningsOn_; // Enable warnings about accessing non-digitized quantities.
272 
273  // Info about tracker module containing stub.
275 
276  // Used to degrade stub bend information.
278 
279  // These module variables are needed only to support the Hybrid stub constructor.
280  // (Otherwise, they could be taken from trackerModule_).
281  unsigned int layerId_;
282  unsigned int layerIdReduced_;
283  float tiltAngle_;
284  float stripPitch_;
286  unsigned int nStrips_;
287  bool psModule_;
288  bool barrel_;
290 
291  const float rejectedStubBend_ = 99999.; // Bend set to this if stub rejected.
292 
293  const float invRoot12 = sqrt(1. / 12.);
294  };
295 
296 } // namespace tmtt
297 #endif
Stub(const Settings *settings, unsigned int idStub, double phi, double r, double z, double bend, unsigned int iphi, double alpha, unsigned int layerId, unsigned int iPhiSec, bool psModule, bool barrel, bool tiltedBarrel, float stripPitch, float stripLength, unsigned int nStrips)
Definition: Stub.cc:21
float stripLength_
Definition: Stub.h:285
float bendCutInFrontend() const
Definition: Stub.h:124
float dphiOverBend_
Definition: Stub.h:243
DigiStage lastDigiStep_
Definition: Stub.h:270
float alpha() const
Definition: Stub.h:119
unsigned int min_qOverPt_bin_
Definition: Stub.h:244
double bendCutExtra() const
Definition: Settings.h:74
bool stubFailedDegradeWindow() const
Definition: Stub.h:159
float sigmaPerp() const
Definition: Stub.h:189
std::array< float, 2 > localU_cluster() const
Definition: Stub.h:153
float stripPitch() const
Definition: Stub.h:205
unsigned int nStrips() const
Definition: Stub.h:209
bool tiltedBarrel() const
Definition: Stub.h:203
std::unique_ptr< DigitalStub > digitalStub_
Definition: Stub.h:269
const Settings * settings_
Definition: Stub.h:233
float eta() const
Definition: Stub.h:111
float sigmaPar() const
Definition: Stub.h:191
unsigned int index_in_vStubs_
Definition: Stub.h:235
unsigned int iphi_
Definition: Stub.h:251
edmNew::DetSet< TTStub< Ref_Phase2TrackerDigi_ > > TTStubDetSet
Definition: Stub.h:34
unsigned int nStrips_
Definition: Stub.h:286
unsigned int layerId_
Definition: Stub.h:281
float sigmaR() const
Definition: Stub.h:186
const TTStubRef & ttStubRef() const
Definition: Stub.h:75
float tiltAngle() const
Definition: Stub.h:184
void setTrackerModule(const TrackerGeometry *trackerGeometry, const TrackerTopology *trackerTopology, const DetId &detId)
std::array< const TP *, 2 > assocTPofCluster_
Definition: Stub.h:267
const DigitalStub * digitalStub() const
Definition: Stub.h:96
std::array< float, 2 > localV_cluster() const
Definition: Stub.h:154
const TrackerModule * trackerModule() const
Definition: Stub.h:78
bool operator==(const Stub &stubOther)
Definition: Stub.h:72
float bend() const
Definition: Stub.h:126
unsigned int iphi() const
Definition: Stub.h:115
Stores association of Truth Particles (TP) to L1 Track-Trigger Clusters.
bool barrel() const
Definition: Stub.h:201
double approxB()
Definition: Stub.cc:304
edm::Ref< TTStubDetSetVec, TTStub< Ref_Phase2TrackerDigi_ > > TTStubRef
unsigned int layerIdReduced() const
Definition: Stub.h:200
Definition: TP.h:23
const DegradeBend * degradeBend_
Definition: Stub.h:277
float phi() const
Definition: Stub.h:107
float trkPhiAtR(float rad) const
Definition: Stub.h:140
unsigned int numMergedBend_
Definition: Stub.h:261
float r_
Definition: Stub.h:240
std::array< float, 2 > localV_cluster_
Definition: Stub.h:249
edmNew::DetSetVector< TTStub< Ref_Phase2TrackerDigi_ > > TTStubDetSetVec
TTStubAssociationMap< Ref_Phase2TrackerDigi_ > TTStubAssMap
Definition: Stub.h:38
T sqrt(T t)
Definition: SSEVec.h:19
unsigned int max_qOverPt_bin_
Definition: Stub.h:245
bool psModule_
Definition: Stub.h:287
void calcQoverPtrange()
Definition: Stub.cc:163
void setDigitizeWarningsOn(bool newVal)
Definition: Stub.h:93
unsigned int layerId() const
Definition: Stub.h:198
double bendCut() const
Definition: Settings.h:70
bool barrel_
Definition: Stub.h:288
float alpha_
Definition: Stub.h:252
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void degradeResolution(float bend, float &degradedBend, unsigned int &num) const
Definition: Stub.cc:248
float tiltAngle_
Definition: Stub.h:283
DigiStage
Definition: Stub.h:89
const TP * assocTP() const
Definition: Stub.h:168
std::set< const TP * > assocTPs_
Definition: Stub.h:265
float z_
Definition: Stub.h:241
float r() const
Definition: Stub.h:108
float qOverPt() const
Definition: Stub.h:132
std::array< bool, 2 > genuineCluster() const
Definition: Stub.h:173
bool genuine() const
Definition: Stub.h:167
double invPtToDphi() const
Definition: Settings.h:397
void setFrontend(const StubKiller *stubKiller)
Definition: Stub.cc:267
std::array< float, 2 > localU_cluster_
Definition: Stub.h:248
float bendCut() const
Definition: Stub.h:128
bool digitizeWarningsOn_
Definition: Stub.h:271
std::array< const TP *, 2 > assocTPofCluster() const
Definition: Stub.h:176
Definition: DetId.h:17
const float rejectedStubBend_
Definition: Stub.h:291
float theta() const
Definition: Stub.h:110
float qOverPtOverBend() const
Definition: Stub.h:148
float sigmaZ() const
Definition: Stub.h:187
float dphiOverBend() const
Definition: Stub.h:146
float qOverPtcut() const
Definition: Stub.h:133
float phi_
Definition: Stub.h:239
unsigned int min_qOverPt_bin() const
Definition: Stub.h:135
float z() const
Definition: Stub.h:109
NOTE: this is needed even if it seems not.
Definition: TTCluster.h:27
TTClusterAssociationMap< Ref_Phase2TrackerDigi_ > TTClusterAssMap
Definition: Stub.h:39
const TP * assocTP_
Definition: Stub.h:264
TTStubRef ttStubRef_
Definition: Stub.h:231
bool frontendPass() const
Definition: Stub.h:157
float bendInFrontend() const
Definition: Stub.h:123
void calcDphiOverBend()
Definition: Stub.cc:314
unsigned int max_qOverPt_bin() const
Definition: Stub.h:136
float stripPitch_
Definition: Stub.h:284
bool psModule() const
Definition: Stub.h:196
edm::Ref< edmNew::DetSetVector< TTCluster< Ref_Phase2TrackerDigi_ > >, TTCluster< Ref_Phase2TrackerDigi_ > > TTClusterRef
Definition: Stub.h:37
unsigned int layerIdReduced_
Definition: Stub.h:282
void fillTruth(const std::map< edm::Ptr< TrackingParticle >, const TP * > &translateTP, const edm::Handle< TTStubAssMap > &mcTruthTTStubHandle, const edm::Handle< TTClusterAssMap > &mcTruthTTClusterHandle)
Definition: Stub.cc:328
float numMergedBend() const
Definition: Stub.h:130
const std::set< const TP * > & assocTPs() const
Definition: Stub.h:164
float bend_
Definition: Stub.h:242
float trkPhiAtRcut(float rad) const
Definition: Stub.h:142
bool frontendPass_
Definition: Stub.h:255
const TrackerModule * trackerModule_
Definition: Stub.h:274
float stripLength() const
Definition: Stub.h:207
void digitize(unsigned int iPhiSec, DigiStage digiStep)
Definition: Stub.cc:203
bool tiltedBarrel_
Definition: Stub.h:289
float phiDiff(float rad, float Pt) const
Definition: Stub.h:138
bool stubFailedDegradeWindow_
Definition: Stub.h:257
Stores association of Truth Particles (TP) to L1 Track-Trigger Stubs.
unsigned int index() const
Definition: Stub.h:101
const float invRoot12
Definition: Stub.h:293
float bendInFrontend_
Definition: Stub.h:259