CMS 3D CMS Logo

Stub.cc
Go to the documentation of this file.
4 
8 
9 #include <cmath>
10 #include <bitset>
11 
12 using namespace std;
13 using namespace trklet;
14 
15 Stub::Stub(Settings const& settings) : settings_(settings) {}
16 
17 Stub::Stub(L1TStub& stub, Settings const& settings, Globals& globals) : settings_(settings) {
18  const string& stubwordhex = stub.stubword();
19 
20  const string stubwordbin = convertHexToBin(stubwordhex);
21 
22  layerdisk_ = stub.layerdisk();
23 
24  int nbendbits = stub.isPSmodule() ? N_BENDBITS_PS : N_BENDBITS_2S;
25 
26  int nalphabits = 0;
27 
28  int nrbits = settings_.nrbitsstub(layerdisk_);
29  int nzbits = settings_.nzbitsstub(layerdisk_);
30  int nphibits = settings_.nphibitsstub(layerdisk_);
31 
32  if (layerdisk_ >= N_LAYER && !stub.isPSmodule()) {
33  nalphabits = settings.nbitsalpha();
34  nrbits = 7;
35  }
36 
37  assert(nbendbits + nalphabits + nrbits + nzbits + nphibits == 36);
38 
39  bitset<32> rbits(stubwordbin.substr(0, nrbits));
40  bitset<32> zbits(stubwordbin.substr(nrbits, nzbits));
41  bitset<32> phibits(stubwordbin.substr(nrbits + nzbits, nphibits));
42  bitset<32> alphabits(stubwordbin.substr(nphibits + nzbits + nrbits, nalphabits));
43  bitset<32> bendbits(stubwordbin.substr(nphibits + nzbits + nrbits + nalphabits, nbendbits));
44 
45  int newbend = bendbits.to_ulong();
46 
47  int newr = rbits.to_ulong();
48  if (layerdisk_ < N_LAYER) {
49  if (newr >= (1 << (nrbits - 1)))
50  newr = newr - (1 << nrbits);
51  }
52 
53  int newz = zbits.to_ulong();
54  if (newz >= (1 << (nzbits - 1)))
55  newz = newz - (1 << nzbits);
56 
57  int newphi = phibits.to_ulong();
58 
59  int newalpha = alphabits.to_ulong();
60  if (newalpha >= (1 << (nalphabits - 1)))
61  newalpha = newalpha - (1 << nalphabits);
62 
63  l1tstub_ = &stub;
64 
65  bend_.set(newbend, nbendbits, true, __LINE__, __FILE__);
66 
67  phi_.set(newphi, nphibits, true, __LINE__, __FILE__);
68  phicorr_.set(newphi, nphibits, true, __LINE__, __FILE__);
69  bool pos = false;
70  if (layerdisk_ >= N_LAYER) {
71  pos = true;
72  int disk = layerdisk_ - N_LAYER + 1;
73  if (stub.z() < 0.0)
74  disk = -disk;
75  disk_.set(disk, 4, false, __LINE__, __FILE__);
76  if (!stub.isPSmodule()) {
77  alpha_.set(newalpha, nalphabits, false, __LINE__, __FILE__);
78  nrbits = 4;
79  }
80  int negdisk = (disk < 0) ? 1 : 0;
81  negdisk_.set(negdisk, 1, true, __LINE__, __FILE__);
82  } else {
83  disk_.set(0, 4, false, __LINE__, __FILE__);
84  layer_.set(layerdisk_, 3, true, __LINE__, __FILE__);
85  }
86  r_.set(newr, nrbits, pos, __LINE__, __FILE__);
87  z_.set(newz, nzbits, false, __LINE__, __FILE__);
88 
89  if (settings.writeMonitorData("StubBend")) {
90  unsigned int nsimtrks = globals.event()->nsimtracks();
91 
92  for (unsigned int isimtrk = 0; isimtrk < nsimtrks; isimtrk++) {
93  const L1SimTrack& simtrk = globals.event()->simtrack(isimtrk);
94  if (stub.tpmatch2(simtrk.trackid())) {
95  double dr = 0.18;
96  double rinv = simtrk.charge() * 0.01 * settings_.c() * settings_.bfield() / simtrk.pt();
97  double pitch = settings_.stripPitch(stub.isPSmodule());
98  double bend = stub.r() * dr * 0.5 * rinv / pitch;
99 
100  globals.ofstream("stubbend.dat") << layerdisk_ << " " << stub.isPSmodule() << " "
101  << simtrk.pt() * simtrk.charge() << " " << bend << " " << newbend << " "
102  << settings.benddecode(newbend, layerdisk_, stub.isPSmodule()) << " "
103  << settings.bendcut(newbend, layerdisk_, stub.isPSmodule()) << endl;
104  }
105  }
106  }
107 }
108 
109 FPGAWord Stub::iphivmFineBins(int VMbits, int finebits) const {
110  unsigned int finephi = (phicorr_.value() >> (phicorr_.nbits() - VMbits - finebits)) & ((1 << finebits) - 1);
111  return FPGAWord(finephi, finebits, true, __LINE__, __FILE__);
112 }
113 
114 unsigned int Stub::phiregionaddress() const {
116  return (iphi << 7) + stubindex_.value();
117 }
118 
121  FPGAWord phiregion(iphi, 3, true, __LINE__, __FILE__);
122  return phiregion.str() + stubindex_.str();
123 }
124 
127  FPGAWord phiregion(iphi, 3, true, __LINE__, __FILE__);
128  return phiregion.str();
129 }
130 
131 void Stub::setAllStubIndex(int nstub) {
132  if (nstub >= (1 << N_BITSMEMADDRESS)) {
133  if (settings_.debugTracklet())
134  edm::LogPrint("Tracklet") << "Warning too large stubindex!";
135  nstub = (1 << N_BITSMEMADDRESS) - 1;
136  }
137 
139 }
140 
141 void Stub::setPhiCorr(int phiCorr) {
142  int iphicorr = phi_.value() - phiCorr;
143 
144  if (iphicorr < 0)
145  iphicorr = 0;
146  if (iphicorr >= (1 << phi_.nbits()))
147  iphicorr = (1 << phi_.nbits()) - 1;
148 
149  phicorr_.set(iphicorr, phi_.nbits(), true, __LINE__, __FILE__);
150 }
151 
152 double Stub::rapprox() const {
153  if (disk_.value() == 0) {
154  int lr = 1 << (8 - settings_.nrbitsstub(layer_.value()));
155  return r_.value() * settings_.kr() * lr + settings_.rmean(layer_.value());
156  }
157  if (!l1tstub_->isPSmodule()) {
158  if (abs(disk_.value()) <= 2)
159  return settings_.rDSSinner(r_.value());
160  else
161  return settings_.rDSSouter(r_.value());
162  }
163  return r_.value() * settings_.kr();
164 }
165 
166 double Stub::zapprox() const {
167  if (disk_.value() == 0) {
168  int lz = 1;
169  if (layer_.value() >= 3) {
170  lz = 16;
171  }
172  return z_.value() * settings_.kz() * lz;
173  }
174  int sign = 1;
175  if (disk_.value() < 0)
176  sign = -1;
177  if (sign < 0) {
178  //Should understand why this is needed to get agreement with integer calculations
179  return (z_.value() + 1) * settings_.kz() + sign * settings_.zmean(abs(disk_.value()) - 1);
180  } else {
181  return z_.value() * settings_.kz() + sign * settings_.zmean(abs(disk_.value()) - 1);
182  }
183 }
184 
185 double Stub::phiapprox(double phimin, double) const {
186  int lphi = 1;
187  if (layer_.value() >= 3) {
188  lphi = 8;
189  }
190  return reco::reduceRange(phimin + phi_.value() * settings_.kphi() / lphi);
191 }
192 
193 unsigned int Stub::layerdisk() const {
194  if (layer_.value() == -1)
195  return N_LAYER - 1 + abs(disk_.value());
196  return layer_.value();
197 }
FPGAWord phi_
Definition: Stub.h:97
double kz() const
Definition: Settings.h:342
double zapprox() const
Definition: Stub.cc:166
FPGAWord alpha_
Definition: Stub.h:98
unsigned int nrbitsstub(unsigned int layerdisk) const
Definition: Settings.h:93
double rDSSinner(unsigned int iBin) const
Definition: Settings.h:183
const FPGAWord & negdisk() const
Definition: Stub.h:67
constexpr T reduceRange(T x)
Definition: deltaPhi.h:18
double phiapprox(double phimin, double) const
Definition: Stub.cc:185
Stub(Settings const &settings)
Definition: Stub.cc:15
FPGAWord r_
Definition: Stub.h:94
FPGAWord layer_
Definition: Stub.h:92
constexpr unsigned int N_BENDBITS_2S
Definition: Settings.h:34
double z() const
Definition: L1TStub.h:59
const FPGAWord & bend() const
Definition: Stub.h:63
std::string phiregionaddressstr() const
Definition: Stub.cc:119
FPGAWord phicorr_
Definition: Stub.h:102
assert(be >=bs)
constexpr unsigned int N_BENDBITS_PS
Definition: Settings.h:33
FPGAWord z_
Definition: Stub.h:95
unsigned int isPSmodule() const
Definition: L1TStub.h:103
const FPGAWord & disk() const
Definition: Stub.h:74
FPGAWord iphivmFineBins(int VMbits, int finebits) const
Definition: Stub.cc:109
FPGAWord disk_
Definition: Stub.h:93
const L1SimTrack & simtrack(int i) const
Definition: SLHCEvent.h:69
unsigned int nbitsallstubs(unsigned int layerdisk) const
Definition: Settings.h:115
unsigned int nzbitsstub(unsigned int layerdisk) const
Definition: Settings.h:91
unsigned int nphibitsstub(unsigned int layerdisk) const
Definition: Settings.h:92
double pt() const
Definition: L1SimTrack.h:26
double rmean(unsigned int iLayer) const
Definition: Settings.h:173
double rDSSouter(unsigned int iBin) const
Definition: Settings.h:186
void setAllStubIndex(int nstub)
Definition: Stub.cc:131
FPGAWord bend_
Definition: Stub.h:100
double bfield() const
Definition: Settings.h:277
double stripPitch(bool isPSmodule) const
Definition: Settings.h:284
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned int phiregionaddress() const
Definition: Stub.cc:114
int value() const
Definition: FPGAWord.h:24
int layerdisk() const
Definition: L1TStub.h:119
bool tpmatch2(int tp) const
Definition: L1TStub.cc:117
FPGAWord negdisk_
Definition: Stub.h:96
bool writeMonitorData(std::string module) const
Definition: Settings.h:118
unsigned int layerdisk() const
Definition: Stub.cc:193
double rapprox() const
Definition: Stub.cc:152
double zmean(unsigned int iDisk) const
Definition: Settings.h:176
int charge() const
Definition: L1SimTrack.h:34
bool debugTracklet() const
Definition: Settings.h:194
int trackid() const
Definition: L1SimTrack.h:24
void setPhiCorr(int phiCorr)
Definition: Stub.cc:141
double kr() const
Definition: Settings.h:344
int nbits() const
Definition: FPGAWord.h:25
double bendcut(unsigned int ibend, unsigned int layerdisk, bool isPSmodule) const
Definition: Settings.h:451
double rinv(double phi1, double phi2, double r1, double r2)
Definition: Util.h:66
const std::string & stubword() const
Definition: L1TStub.h:123
void set(int value, int nbits, bool positive=true, int line=-1, const char *file=nullptr)
Definition: FPGAWord.cc:14
SLHCEvent *& event()
Definition: Globals.h:36
unsigned int nsimtracks() const
Definition: SLHCEvent.h:67
L1TStub * l1tstub_
Definition: Stub.h:106
double c() const
Definition: Settings.h:224
double r() const
Definition: L1TStub.h:60
std::ofstream & ofstream(std::string fname)
Definition: Globals.cc:44
Settings const & settings_
Definition: Stub.h:107
int nbitsalpha() const
Definition: Settings.h:229
double kphi() const
Definition: Settings.h:338
std::string convertHexToBin(const std::string &stubwordhex)
Definition: Util.h:75
std::string str() const
Definition: FPGAWord.cc:54
std::string phiregionstr() const
Definition: Stub.cc:125
FPGAWord stubindex_
Definition: Stub.h:104
constexpr unsigned int N_BITSMEMADDRESS
Definition: Settings.h:43
double benddecode(unsigned int ibend, unsigned int layerdisk, bool isPSmodule) const
Definition: Settings.h:443
unsigned int layerdisk_
Definition: Stub.h:90
constexpr int N_LAYER
Definition: Settings.h:25