File indexing completed on 2023-03-17 11:13:41
0001 #include "L1Trigger/TrackFindingTMTT/interface/HTrphi.h"
0002 #include "L1Trigger/TrackFindingTMTT/interface/InputData.h"
0003 #include "L1Trigger/TrackFindingTMTT/interface/TP.h"
0004 #include "L1Trigger/TrackFindingTMTT/interface/L1fittedTrack.h"
0005 #include "L1Trigger/TrackFindingTMTT/interface/Settings.h"
0006 #include "L1Trigger/TrackFindingTMTT/interface/PrintL1trk.h"
0007 #include "FWCore/ServiceRegistry/interface/Service.h"
0008
0009 #include "DataFormats/Math/interface/deltaPhi.h"
0010 #include "FWCore/Utilities/interface/Exception.h"
0011
0012 #include <vector>
0013 #include <limits>
0014 #include <atomic>
0015 #include <sstream>
0016 #include <mutex>
0017
0018 using namespace std;
0019
0020 namespace tmtt {
0021
0022
0023
0024
0025
0026
0027
0028
0029 HTrphi::HTrphi(const Settings* settings,
0030 unsigned int iPhiSec,
0031 unsigned int iEtaReg,
0032 float etaMinSector,
0033 float etaMaxSector,
0034 float phiCentreSector,
0035 HTrphi::ErrorMonitor* errMon)
0036 : HTbase(settings, iPhiSec, iEtaReg, settings->houghNbinsPt(), settings->houghNbinsPhi()),
0037 invPtToDphi_((settings->invPtToDphi())),
0038 shape_(static_cast<HTshape>(settings->shape())),
0039
0040
0041 maxAbsQoverPtAxis_(1. / settings->houghMinPt()),
0042 nBinsQoverPtAxis_(settings->houghNbinsPt()),
0043 binSizeQoverPtAxis_(2 * maxAbsQoverPtAxis_ / nBinsQoverPtAxis_),
0044
0045
0046
0047 chosenRofPhi_(settings->chosenRofPhi()),
0048 phiCentreSector_(phiCentreSector),
0049 maxAbsPhiTrkAxis_(M_PI / float(settings->numPhiSectors())),
0050 nBinsPhiTrkAxis_(settings->houghNbinsPhi()),
0051 binSizePhiTrkAxis_(2 * maxAbsPhiTrkAxis_ / nBinsPhiTrkAxis_),
0052 errMon_(errMon) {
0053
0054 if (shape_ != HTshape::square)
0055 binSizeQoverPtAxis_ = 2. * maxAbsQoverPtAxis_ / (nBinsQoverPtAxis_ - 1.);
0056 if (shape_ == HTshape::hexagon)
0057 binSizePhiTrkAxis_ = 2. * maxAbsPhiTrkAxis_ / (nBinsPhiTrkAxis_ - 1. / 6.);
0058 else if (shape_ == HTshape::diamond)
0059 binSizePhiTrkAxis_ = 2. * maxAbsPhiTrkAxis_ / (nBinsPhiTrkAxis_ - 1. / 2.);
0060
0061
0062
0063
0064 enableMerge2x2_ = (settings->enableMerge2x2() || settings->miniHTstage());
0065 if (settings->miniHTstage()) {
0066
0067 minInvPtToMerge2x2_ = 0.;
0068 } else {
0069
0070 minInvPtToMerge2x2_ = 1. / (settings->maxPtToMerge2x2());
0071 if (minInvPtToMerge2x2_ > maxAbsQoverPtAxis_)
0072 enableMerge2x2_ = false;
0073 }
0074
0075
0076
0077
0078
0079 if (enableMerge2x2_ && (nBinsQoverPtAxis_ % 2 != 0 || nBinsPhiTrkAxis_ % 2 != 0))
0080 throw cms::Exception("BadConfig") << "HTrphi: You are not allowed to set EnableMerge2x2 or MiniHTstage = True if "
0081 "you have an odd number of bins "
0082 "in r-phi HT array "
0083 << nBinsQoverPtAxis_ << " " << nBinsPhiTrkAxis_;
0084
0085
0086
0087
0088 killSomeHTCellsRphi_ = settings->killSomeHTCellsRphi();
0089
0090
0091 nReceivedStubs_ = 0;
0092 busyInputSectorKill_ = settings_->busyInputSectorKill();
0093 busySectorKill_ = settings_->busySectorKill();
0094
0095 busyInputSectorNumStubs_ = settings_->busyInputSectorNumStubs();
0096
0097 busySectorNumStubs_ = settings_->busySectorNumStubs();
0098
0099 busySectorMbinRanges_ = settings_->busySectorMbinRanges();
0100
0101 busySectorMbinOrder_ = settings_->busySectorMbinOrder();
0102
0103 busySectorUseMbinRanges_ = (not busySectorMbinRanges_.empty());
0104 busySectorUseMbinOrder_ = (not busySectorMbinOrder_.empty());
0105
0106 bool rescaleMbins = false;
0107 if (busySectorUseMbinRanges_) {
0108
0109
0110
0111 unsigned int nTotalBins = 0;
0112 for (unsigned int j = 0; j < busySectorMbinRanges_.size(); j++) {
0113 nTotalBins += busySectorMbinRanges_[j];
0114 }
0115 rescaleMbins = (nTotalBins != nBinsQoverPtAxis_);
0116
0117 if (rescaleMbins && busySectorUseMbinOrder_)
0118 throw cms::Exception("BadConfig") << "HTrphi: BusySectorUserMbin error";
0119 float rescaleFactor = rescaleMbins ? float(nBinsQoverPtAxis_) / float(nTotalBins) : 1.;
0120
0121 busySectorMbinLow_.resize(busySectorMbinRanges_.size());
0122 busySectorMbinHigh_.resize(busySectorMbinRanges_.size());
0123 float mBinSum = 0.;
0124 for (unsigned int i = 0; i < busySectorMbinRanges_.size(); i++) {
0125 busySectorMbinLow_[i] = std::round(mBinSum);
0126 busySectorMbinHigh_[i] = std::round(mBinSum + rescaleFactor * busySectorMbinRanges_[i]) - 1;
0127 mBinSum += rescaleFactor * busySectorMbinRanges_[i];
0128 }
0129 }
0130
0131 for (unsigned int i = 0; i < nBinsQoverPtAxis_; i++) {
0132 for (unsigned int j = 0; j < nBinsPhiTrkAxis_; j++) {
0133 pair<float, float> helix = this->helix2Dconventional(i, j);
0134 float qOverPt = helix.first;
0135
0136 bool mergedCell = false;
0137 if (enableMerge2x2_ && this->mergedCell(i, j))
0138 mergedCell = true;
0139
0140 HTbase::htArray_(i, j) =
0141 std::make_unique<HTcell>(settings, iPhiSec, iEtaReg, etaMinSector, etaMaxSector, qOverPt, i, mergedCell);
0142 }
0143 }
0144
0145 std::stringstream text;
0146 text << "\n";
0147 text << "=== R-PHI HOUGH TRANSFORM AXES RANGES: abs(q/Pt) < " << maxAbsQoverPtAxis_ << " & abs(track-phi) < "
0148 << maxAbsPhiTrkAxis_ << " ===\n";
0149 text << "=== R-PHI HOUGH TRANSFORM ARRAY SIZE: q/Pt bins = " << nBinsQoverPtAxis_
0150 << " & track-phi bins = " << nBinsPhiTrkAxis_ << " ===\n";
0151 text << "=== R-PHI HOUGH TRANSFORM BIN SIZE: BIN(q/Pt) = " << binSizeQoverPtAxis_
0152 << " & BIN(track-phi) = " << binSizePhiTrkAxis_ << " ===\n\n";
0153 if (busySectorKill_ && busySectorUseMbinRanges_ && rescaleMbins) {
0154 text << "=== R-PHI HOUGH TRANSFORM WARNING: Rescaled m bin ranges specified by cfg parameter "
0155 "BusySectorMbinRanges, as they were inconsistent with total number of m bins in HT.\n";
0156 text << "=== Rescaled values for BusySectorMbinRanges =";
0157 for (unsigned int i = 0; i < busySectorMbinRanges_.size(); i++) {
0158 text << " " << (busySectorMbinHigh_[i] - busySectorMbinLow_[i] + 1);
0159 }
0160 }
0161 text << "\n";
0162 static std::once_flag printOnce;
0163 std::call_once(
0164 printOnce, [](string t) { PrintL1trk() << t; }, text.str());
0165
0166
0167 cellCenters_.clear();
0168 for (unsigned int m = 0; m < nBinsQoverPtAxis_; m++) {
0169 std::vector<std::pair<float, float> > binCenters;
0170 for (unsigned int c = 0; c < nBinsPhiTrkAxis_; c++)
0171 binCenters.push_back(this->helix2Dhough(m, c));
0172 cellCenters_.push_back(binCenters);
0173 }
0174 }
0175
0176
0177
0178
0179 void HTrphi::store(Stub* stub, const vector<bool>& inEtaSubSecs) {
0180
0181 if ((!busyInputSectorKill_) || (nReceivedStubs_ < busyInputSectorNumStubs_)) {
0182 nReceivedStubs_++;
0183
0184 unsigned int jPhiTrkBinMinLast = 0;
0185 unsigned int jPhiTrkBinMaxLast = 99999;
0186
0187
0188 for (unsigned int i = 0; i < nBinsQoverPtAxis_; i++) {
0189 if (shape_ == HTshape::square) {
0190
0191
0192
0193 pair<unsigned int, unsigned int> jRange = this->iPhiRange(stub, i);
0194 unsigned int jPhiTrkBinMin = jRange.first;
0195 unsigned int jPhiTrkBinMax = jRange.second;
0196
0197
0198 for (unsigned int j = jPhiTrkBinMin; j <= jPhiTrkBinMax; j++) {
0199 bool canStoreStub = true;
0200 unsigned int iStore = i;
0201 unsigned int jStore = j;
0202
0203
0204
0205 if (enableMerge2x2_) {
0206
0207 if (this->mergedCell(i, j)) {
0208
0209
0210 if (i % 2 == 1)
0211 iStore = i - 1;
0212 if (j % 2 == 1)
0213 jStore = j - 1;
0214
0215 if (HTbase::htArray_(iStore, jStore)->stubStoredInCell(stub))
0216 canStoreStub = false;
0217 }
0218 }
0219
0220 if (canStoreStub)
0221 HTbase::htArray_(iStore, jStore)->store(stub, inEtaSubSecs);
0222 }
0223
0224
0225 if (errMon_ != nullptr) {
0226 this->countFirmwareErrors(i, jPhiTrkBinMin, jPhiTrkBinMax, jPhiTrkBinMinLast, jPhiTrkBinMaxLast);
0227 jPhiTrkBinMinLast = jPhiTrkBinMin;
0228 jPhiTrkBinMaxLast = jPhiTrkBinMax;
0229 }
0230
0231 } else {
0232
0233
0234 if (shape_ == HTshape::diamond) {
0235
0236
0237 float qOverPtBin = -maxAbsQoverPtAxis_ + i * binSizeQoverPtAxis_;
0238 float phiTrk = reco::deltaPhi(stub->phi(), phiCentreSector_) +
0239 invPtToDphi_ * qOverPtBin * (stub->r() - chosenRofPhi_) + maxAbsPhiTrkAxis_;
0240 if (i % 2 == 0)
0241 phiTrk += binSizePhiTrkAxis_ / 2.;
0242 unsigned int binCenter = std::floor(phiTrk / binSizePhiTrkAxis_);
0243 if (binCenter < nBinsPhiTrkAxis_)
0244 HTbase::htArray_(i, binCenter)->store(stub, inEtaSubSecs);
0245
0246 } else if (shape_ == HTshape::hexagon) {
0247
0248
0249 float qOverPtBin = -maxAbsQoverPtAxis_ + i * binSizeQoverPtAxis_;
0250 float qOverPtBinVar = binSizeQoverPtAxis_;
0251 float phiTrk = reco::deltaPhi(stub->phi(), phiCentreSector_) +
0252 invPtToDphi_ * qOverPtBin * (stub->r() - chosenRofPhi_) + maxAbsPhiTrkAxis_;
0253 float phiTrkVar = invPtToDphi_ * qOverPtBinVar * std::abs(stub->r() - chosenRofPhi_);
0254 float phiTrkMin = phiTrk - phiTrkVar;
0255 float phiTrkMax = phiTrk + phiTrkVar;
0256 if (i % 2 == 0)
0257 phiTrk += binSizePhiTrkAxis_ / 6.;
0258 else {
0259 phiTrk -= binSizePhiTrkAxis_ / 3.;
0260 phiTrkMin -= binSizePhiTrkAxis_ / 2.;
0261 phiTrkMax -= binSizePhiTrkAxis_ / 2.;
0262 }
0263 unsigned int iCenter = std::floor(phiTrk / binSizePhiTrkAxis_ * 3.);
0264 unsigned int iMin = std::floor(phiTrkMin / binSizePhiTrkAxis_ * 3.);
0265 unsigned int iMax = std::floor(phiTrkMax / binSizePhiTrkAxis_ * 3.);
0266 std::pair<bool, unsigned int> binCenter;
0267 std::pair<bool, unsigned int> binMin;
0268 std::pair<bool, unsigned int> binMax;
0269 binCenter.second = iCenter / 3;
0270 binMin.second = iMin / 3;
0271 binMax.second = iMax / 3;
0272 binCenter.first = !(iCenter % 3 == 2);
0273 binMin.first = (iMin % 3 == 0);
0274 binMax.first = (iMax % 3 == 0);
0275 if (binCenter.first && binCenter.second < nBinsPhiTrkAxis_)
0276 HTbase::htArray_(i, binCenter.second)->store(stub, inEtaSubSecs);
0277 else if (binMin.first && binMin.second < nBinsPhiTrkAxis_)
0278 HTbase::htArray_(i, binMin.second)->store(stub, inEtaSubSecs);
0279 else if (binMax.first && binMax.second < nBinsPhiTrkAxis_)
0280 HTbase::htArray_(i, binMax.second)->store(stub, inEtaSubSecs);
0281
0282 } else if (shape_ == HTshape::brick) {
0283
0284
0285 float qOverPtBin = -maxAbsQoverPtAxis_ + i * binSizeQoverPtAxis_;
0286 float qOverPtBinVar = binSizeQoverPtAxis_;
0287 float phiTrk = reco::deltaPhi(stub->phi(), phiCentreSector_) +
0288 invPtToDphi_ * qOverPtBin * (stub->r() - chosenRofPhi_) + maxAbsPhiTrkAxis_;
0289 float phiTrkVar = invPtToDphi_ * qOverPtBinVar * std::abs(stub->r() - chosenRofPhi_);
0290 float phiTrkMin = phiTrk - phiTrkVar;
0291 float phiTrkMax = phiTrk + phiTrkVar;
0292 unsigned int iMin = std::floor(phiTrkMin / binSizePhiTrkAxis_ * 2.);
0293 unsigned int iMax = std::floor(phiTrkMax / binSizePhiTrkAxis_ * 2.);
0294 std::pair<bool, unsigned int> binMin;
0295 std::pair<bool, unsigned int> binMax;
0296 binMin.second = iMin / 2;
0297 binMax.second = iMax / 2;
0298 binMin.first = (iMin % 2 == i % 2);
0299 binMax.first = (iMax % 2 == i % 2);
0300 if (binMin.first && binMin.second < nBinsPhiTrkAxis_)
0301 HTbase::htArray_(i, binMin.second)->store(stub, inEtaSubSecs);
0302 else if (binMax.first && binMax.second < nBinsPhiTrkAxis_)
0303 HTbase::htArray_(i, binMax.second)->store(stub, inEtaSubSecs);
0304 }
0305 }
0306 }
0307
0308
0309 if (errMon_ != nullptr) {
0310 errMon_->maxLineGradient = max(errMon_->maxLineGradient.load(), this->calcLineGradArray(stub->r()));
0311 }
0312 }
0313 }
0314
0315
0316
0317 unsigned int HTrphi::getMbinRange(const L1track2D& trk) const {
0318 if (busySectorUseMbinRanges_) {
0319 unsigned int mBin = trk.cellLocationHT().first;
0320 unsigned int mBinOrder;
0321 if (busySectorUseMbinOrder_) {
0322
0323 mBinOrder = 99999;
0324 for (unsigned int k = 0; k < busySectorMbinOrder_.size(); k++) {
0325 if (mBin == busySectorMbinOrder_[k])
0326 mBinOrder = k;
0327 }
0328 if (mBinOrder == 99999)
0329 throw cms::Exception("LogicError") << "HTrphi::getMbinRange() mBinOrder calculation wrong.";
0330 } else {
0331
0332 mBinOrder = mBin;
0333 }
0334 for (unsigned int i = 0; i < busySectorMbinRanges_.size(); i++) {
0335 if (mBinOrder >= busySectorMbinLow_[i] && mBinOrder <= busySectorMbinHigh_[i])
0336 return i;
0337 }
0338 throw cms::Exception("LogicError") << "HTrphi::getMbinRange() messed up";
0339 } else {
0340 return 0;
0341 }
0342 }
0343
0344
0345
0346
0347
0348 pair<unsigned int, unsigned int> HTrphi::iPhiRange(const Stub* stub, unsigned int iQoverPtBin, bool debug) const {
0349
0350 float qOverPtBin = -maxAbsQoverPtAxis_ + (iQoverPtBin + 0.5) * binSizeQoverPtAxis_;
0351
0352 float qOverPtBinVar = 0.5 * binSizeQoverPtAxis_;
0353
0354
0355
0356
0357
0358 float phiTrk = stub->phi() + invPtToDphi_ * qOverPtBin * (stub->r() - chosenRofPhi_);
0359
0360
0361
0362 float phiTrkVar = invPtToDphi_ * qOverPtBinVar * std::abs(stub->r() - chosenRofPhi_);
0363 float phiTrkMin = phiTrk - phiTrkVar;
0364 float phiTrkMax = phiTrk + phiTrkVar;
0365
0366 float deltaPhiMin = reco::deltaPhi(phiTrkMin, phiCentreSector_);
0367 float deltaPhiMax = reco::deltaPhi(phiTrkMax, phiCentreSector_);
0368 pair<float, float> phiTrkRange(deltaPhiMin, deltaPhiMax);
0369
0370
0371 pair<unsigned int, unsigned int> iPhiTrkBinRange = this->HTbase::convertCoordRangeToBinRange(
0372 phiTrkRange, nBinsPhiTrkAxis_, (-maxAbsPhiTrkAxis_), binSizePhiTrkAxis_, killSomeHTCellsRphi_);
0373
0374 return iPhiTrkBinRange;
0375 }
0376
0377
0378
0379 void HTrphi::countFirmwareErrors(unsigned int iQoverPtBin,
0380 unsigned int jPhiTrkBinMin,
0381 unsigned int jPhiTrkBinMax,
0382 unsigned int jPhiTrkBinMinLast,
0383 unsigned int jPhiTrkBinMaxLast) {
0384
0385 if (jPhiTrkBinMax >= jPhiTrkBinMin) {
0386
0387
0388 bool OK_a = (jPhiTrkBinMin + 1 >= jPhiTrkBinMinLast) && (jPhiTrkBinMax <= jPhiTrkBinMaxLast + 1);
0389
0390 bool OK_b = (jPhiTrkBinMax - jPhiTrkBinMin + 1 <= 2);
0391
0392 if (!OK_a)
0393 errMon_->numErrorsTypeA++;
0394 if (!OK_b)
0395 errMon_->numErrorsTypeB++;
0396 errMon_->numErrorsNorm++;
0397 }
0398 }
0399
0400
0401
0402
0403
0404 pair<float, float> HTrphi::helix2Dhough(unsigned int i, unsigned int j) const {
0405 unsigned int qOverPtBin = i;
0406 unsigned int phiTrkBin = j;
0407
0408
0409 bool merged = false;
0410 if (enableMerge2x2_) {
0411
0412 if (this->mergedCell(i, j)) {
0413 merged = true;
0414
0415
0416 if (i % 2 == 1)
0417 qOverPtBin = i - 1;
0418 if (j % 2 == 1)
0419 phiTrkBin = j - 1;
0420 }
0421 }
0422
0423 float qOverPtBinCenter = .5;
0424 float phiTrkBinCenter = .5;
0425
0426 if (shape_ != HTshape::square) {
0427 qOverPtBinCenter = 0.;
0428
0429 float evenPhiPos = 0., oddPhiPos = 0.;
0430 if (shape_ == HTshape::hexagon) {
0431 evenPhiPos = 1. / 6.;
0432 oddPhiPos = 2. / 3.;
0433 } else if (shape_ == HTshape::diamond) {
0434 evenPhiPos = 0.;
0435 oddPhiPos = 0.5;
0436 } else if (shape_ == HTshape::brick) {
0437 evenPhiPos = 0.25;
0438 oddPhiPos = 0.75;
0439 }
0440 phiTrkBinCenter = (qOverPtBin % 2 == 0) ? evenPhiPos : oddPhiPos;
0441 }
0442
0443 float qOverPt = -maxAbsQoverPtAxis_ + (qOverPtBin + qOverPtBinCenter) * binSizeQoverPtAxis_;
0444 float phiTrk = -maxAbsPhiTrkAxis_ + (phiTrkBin + phiTrkBinCenter) * binSizePhiTrkAxis_;
0445
0446 if (merged) {
0447 qOverPt += 0.5 * binSizeQoverPtAxis_;
0448 phiTrk += 0.5 * binSizePhiTrkAxis_;
0449 }
0450
0451
0452 phiTrk = reco::deltaPhi(phiTrk + phiCentreSector_, 0.);
0453 return pair<float, float>(qOverPt, phiTrk);
0454 }
0455
0456
0457
0458
0459
0460 pair<float, float> HTrphi::helix2Dconventional(unsigned int i, unsigned int j) const {
0461
0462 pair<float, float> helix2Dht = this->helix2Dhough(i, j);
0463
0464 float qOverPt = helix2Dht.first;
0465
0466 float phi0 = reco::deltaPhi(helix2Dht.second + invPtToDphi_ * chosenRofPhi_ * qOverPt, 0.);
0467 return pair<float, float>(qOverPt, phi0);
0468 }
0469
0470
0471
0472
0473 pair<unsigned int, unsigned int> HTrphi::trueCell(const TP* tp) const {
0474
0475 float qOverPt = tp->qOverPt();
0476 float phiTrk = tp->trkPhiAtR(chosenRofPhi_);
0477
0478 float deltaPhi = reco::deltaPhi(phiTrk, phiCentreSector_);
0479
0480 int iQoverPt = floor((qOverPt - (-maxAbsQoverPtAxis_)) / binSizeQoverPtAxis_);
0481 int iPhiTrk = floor((deltaPhi - (-maxAbsPhiTrkAxis_)) / binSizePhiTrkAxis_);
0482
0483 if (iQoverPt >= 0 && iQoverPt < int(nBinsQoverPtAxis_) && iPhiTrk >= 0 && iPhiTrk < int(nBinsPhiTrkAxis_)) {
0484
0485
0486
0487
0488 ;
0489 } else {
0490
0491 if (iQoverPt < 0)
0492 iQoverPt = 0;
0493 if (iQoverPt >= int(nBinsQoverPtAxis_))
0494 iQoverPt = nBinsQoverPtAxis_ - 1;
0495 if (iPhiTrk < 0)
0496 iPhiTrk = 0;
0497 if (iPhiTrk >= int(nBinsPhiTrkAxis_))
0498 iPhiTrk = nBinsPhiTrkAxis_ - 1;
0499 }
0500 return pair<unsigned int, unsigned int>(iQoverPt, iPhiTrk);
0501 }
0502
0503
0504
0505
0506
0507 pair<unsigned int, unsigned int> HTrphi::cell(const L1fittedTrack* fitTrk) const {
0508 bool beamConstraint = fitTrk->done_bcon();
0509
0510 float qOverPt = beamConstraint ? fitTrk->qOverPt_bcon() : fitTrk->qOverPt();
0511
0512 float phiTrk = fitTrk->phiAtChosenR(beamConstraint);
0513
0514 float deltaPhi = reco::deltaPhi(phiTrk, phiCentreSector_);
0515
0516 int iQoverPt = 999999;
0517 int iPhiTrk = 999999;
0518
0519 if (shape_ == HTshape::square) {
0520
0521
0522 iQoverPt = floor((qOverPt - (-maxAbsQoverPtAxis_)) / binSizeQoverPtAxis_);
0523 iPhiTrk = floor((deltaPhi - (-maxAbsPhiTrkAxis_)) / binSizePhiTrkAxis_);
0524
0525
0526 if (iQoverPt >= 0 && iQoverPt < int(nBinsQoverPtAxis_) && iPhiTrk >= 0 && iPhiTrk < int(nBinsPhiTrkAxis_)) {
0527
0528
0529
0530
0531 ;
0532 } else {
0533
0534 if (iQoverPt < 0)
0535 iQoverPt = 0;
0536 if (iQoverPt >= int(nBinsQoverPtAxis_))
0537 iQoverPt = nBinsQoverPtAxis_ - 1;
0538 if (iPhiTrk < 0)
0539 iPhiTrk = 0;
0540 if (iPhiTrk >= int(nBinsPhiTrkAxis_))
0541 iPhiTrk = nBinsPhiTrkAxis_ - 1;
0542 }
0543
0544 } else {
0545
0546
0547 float minD = std::numeric_limits<float>::infinity();
0548 float d(0);
0549 unsigned int m(0);
0550 for (const auto& binCenters : cellCenters_) {
0551 unsigned int c(0);
0552 for (auto cellCenter : binCenters) {
0553 d = std::pow((cellCenter.first - qOverPt) / (float)binSizeQoverPtAxis_, 2) +
0554 std::pow((cellCenter.second - phiTrk) / (float)binSizePhiTrkAxis_, 2);
0555 if (d < minD) {
0556 minD = d;
0557 iQoverPt = m;
0558 iPhiTrk = c;
0559 }
0560 c++;
0561 }
0562 m++;
0563 }
0564
0565 if (iQoverPt < 0)
0566 iQoverPt = 0;
0567 if (iQoverPt >= int(nBinsQoverPtAxis_))
0568 iQoverPt = nBinsQoverPtAxis_ - 1;
0569 if (iPhiTrk < 0)
0570 iPhiTrk = 0;
0571 if (iPhiTrk >= int(nBinsPhiTrkAxis_))
0572 iPhiTrk = nBinsPhiTrkAxis_ - 1;
0573 }
0574
0575 return pair<unsigned int, unsigned int>(iQoverPt, iPhiTrk);
0576 }
0577
0578
0579
0580
0581 bool HTrphi::mergedCell(unsigned int iQoverPtBin, unsigned int jPhiTrkBin) const {
0582 bool merge = false;
0583
0584 if (enableMerge2x2_) {
0585 unsigned int i = iQoverPtBin;
0586
0587
0588
0589 float fMergeBins = (maxAbsQoverPtAxis_ - minInvPtToMerge2x2_) / (2. * binSizeQoverPtAxis_);
0590
0591 unsigned int numQoverPtBinsToMerge = 2 * min((unsigned int)(std::round(fMergeBins)), (nBinsQoverPtAxis_ / 4));
0592 const float small = 0.001;
0593 if (minInvPtToMerge2x2_ < small && (unsigned int)(std::round(2. * fMergeBins)) % 2 == 1)
0594 numQoverPtBinsToMerge++;
0595 unsigned int iB = (nBinsQoverPtAxis_ - 1) - i;
0596 if (min(i, iB) < numQoverPtBinsToMerge)
0597 merge = true;
0598 }
0599
0600 return merge;
0601 }
0602
0603
0604
0605 float HTrphi::calcLineGradArray(float r) const {
0606 float grad = std::abs(invPtToDphi_ * (r - chosenRofPhi_));
0607
0608 grad *= binSizeQoverPtAxis_ / binSizePhiTrkAxis_;
0609 if (shape_ == HTshape::hexagon)
0610 grad *= 3.;
0611 else if (shape_ == HTshape::diamond)
0612 grad *= 2.;
0613 else if (shape_ == HTshape::brick)
0614 grad *= 4.;
0615 return grad;
0616 }
0617
0618
0619
0620
0621 list<L1track2D> HTrphi::killTracksBusySec(const list<L1track2D>& tracks) const {
0622 list<L1track2D> outTracks;
0623
0624 if (busySectorKill_) {
0625 unsigned int nStubsOut = 0;
0626
0627 vector<unsigned int> nStubsOutInRange(busySectorMbinRanges_.size(), 0);
0628
0629 for (const L1track2D& trk : tracks) {
0630 bool keep = true;
0631 unsigned int nStubs = trk.numStubs();
0632 if (busySectorUseMbinRanges_) {
0633 unsigned int mBinRange = this->getMbinRange(trk);
0634 nStubsOutInRange[mBinRange] += nStubs;
0635 if (nStubsOutInRange[mBinRange] > busySectorNumStubs_)
0636 keep = false;
0637 } else {
0638 nStubsOut += nStubs;
0639 if (nStubsOut > busySectorNumStubs_)
0640 keep = false;
0641 }
0642
0643 if (keep)
0644 outTracks.push_back(trk);
0645 }
0646
0647 } else {
0648 outTracks = tracks;
0649 }
0650
0651 return outTracks;
0652 }
0653
0654
0655
0656
0657
0658 vector<unsigned int> HTrphi::rowOrder(unsigned int numRows) const {
0659 vector<unsigned int> iOrder;
0660
0661
0662 const bool oddNumRows = (numRows % 2 == 1);
0663
0664
0665 if (oddNumRows) {
0666 unsigned int middleRow = (numRows - 1) / 2;
0667 iOrder.push_back(middleRow);
0668 for (unsigned int i = 1; i <= (numRows - 1) / 2; i++) {
0669 iOrder.push_back(middleRow - i);
0670 iOrder.push_back(middleRow + i);
0671 }
0672 } else {
0673 unsigned int startRowPos = numRows / 2;
0674 unsigned int startRowNeg = startRowPos - 1;
0675 for (unsigned int i = 0; i < numRows / 2; i++) {
0676 iOrder.push_back(startRowNeg - i);
0677 iOrder.push_back(startRowPos + i);
0678 }
0679 }
0680
0681 return iOrder;
0682 }
0683
0684 }