Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:22:03

0001 #include "L1Trigger/TrackFindingTracklet/interface/TrackletCalculatorBase.h"
0002 #include "L1Trigger/TrackFindingTracklet/interface/Tracklet.h"
0003 #include "L1Trigger/TrackFindingTracklet/interface/Stub.h"
0004 #include "L1Trigger/TrackFindingTracklet/interface/Globals.h"
0005 #include "L1Trigger/TrackFindingTracklet/interface/HistBase.h"
0006 #include "L1Trigger/TrackFindingTracklet/interface/IMATH_TrackletCalculator.h"
0007 #include "L1Trigger/TrackFindingTracklet/interface/IMATH_TrackletCalculatorDisk.h"
0008 #include "L1Trigger/TrackFindingTracklet/interface/IMATH_TrackletCalculatorOverlap.h"
0009 
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include "FWCore/Utilities/interface/Exception.h"
0012 #include "DataFormats/Math/interface/deltaPhi.h"
0013 #include "L1Trigger/L1TCommon/interface/BitShift.h"
0014 
0015 using namespace std;
0016 using namespace trklet;
0017 
0018 TrackletCalculatorBase::TrackletCalculatorBase(string name, Settings const& settings, Globals* global)
0019     : ProcessBase(name, settings, global) {}
0020 
0021 void TrackletCalculatorBase::exacttracklet(double r1,
0022                                            double z1,
0023                                            double phi1,
0024                                            double r2,
0025                                            double z2,
0026                                            double phi2,
0027                                            double,
0028                                            double& rinv,
0029                                            double& phi0,
0030                                            double& t,
0031                                            double& z0,
0032                                            double phiproj[N_LAYER - 2],
0033                                            double zproj[N_LAYER - 2],
0034                                            double phider[N_LAYER - 2],
0035                                            double zder[N_LAYER - 2],
0036                                            double phiprojdisk[N_DISK],
0037                                            double rprojdisk[N_DISK],
0038                                            double phiderdisk[N_DISK],
0039                                            double rderdisk[N_DISK]) {
0040   double deltaphi = reco::reduceRange(phi1 - phi2);
0041 
0042   double dist = sqrt(r2 * r2 + r1 * r1 - 2 * r1 * r2 * cos(deltaphi));
0043 
0044   rinv = 2 * sin(deltaphi) / dist;
0045 
0046   double phi1tmp = phi1 - phimin_;
0047 
0048   phi0 = reco::reduceRange(phi1tmp + asin(0.5 * r1 * rinv));
0049 
0050   double rhopsi1 = 2 * asin(0.5 * r1 * rinv) / rinv;
0051   double rhopsi2 = 2 * asin(0.5 * r2 * rinv) / rinv;
0052 
0053   t = (z1 - z2) / (rhopsi1 - rhopsi2);
0054 
0055   z0 = z1 - t * rhopsi1;
0056 
0057   for (unsigned int i = 0; i < N_LAYER - 2; i++) {
0058     exactproj(settings_.rmean(settings_.projlayers(iSeed_, i) - 1),
0059               rinv,
0060               phi0,
0061               t,
0062               z0,
0063               phiproj[i],
0064               zproj[i],
0065               phider[i],
0066               zder[i]);
0067   }
0068 
0069   for (unsigned int i = 0; i < N_DISK; i++) {
0070     exactprojdisk(settings_.zmean(i), rinv, phi0, t, z0, phiprojdisk[i], rprojdisk[i], phiderdisk[i], rderdisk[i]);
0071   }
0072 }
0073 
0074 void TrackletCalculatorBase::exacttrackletdisk(double r1,
0075                                                double z1,
0076                                                double phi1,
0077                                                double r2,
0078                                                double z2,
0079                                                double phi2,
0080                                                double,
0081                                                double& rinv,
0082                                                double& phi0,
0083                                                double& t,
0084                                                double& z0,
0085                                                double phiprojLayer[N_PSLAYER],  //=3 (project to PS barrel layers only)
0086                                                double zprojLayer[N_PSLAYER],
0087                                                double phiderLayer[N_PSLAYER],
0088                                                double zderLayer[N_PSLAYER],
0089                                                double phiproj[N_DISK - 2],  //=3 (max project to 3 other disks)
0090                                                double rproj[N_DISK - 2],
0091                                                double phider[N_DISK - 2],
0092                                                double rder[N_DISK - 2]) {
0093   double deltaphi = reco::reduceRange(phi1 - phi2);
0094 
0095   double dist = sqrt(r2 * r2 + r1 * r1 - 2 * r1 * r2 * cos(deltaphi));
0096 
0097   rinv = 2 * sin(deltaphi) / dist;
0098 
0099   double phi1tmp = phi1 - phimin_;
0100 
0101   phi0 = reco::reduceRange(phi1tmp + asin(0.5 * r1 * rinv));
0102 
0103   double rhopsi1 = 2 * asin(0.5 * r1 * rinv) / rinv;
0104   double rhopsi2 = 2 * asin(0.5 * r2 * rinv) / rinv;
0105 
0106   t = (z1 - z2) / (rhopsi1 - rhopsi2);
0107 
0108   z0 = z1 - t * rhopsi1;
0109 
0110   for (unsigned int i = 0; i < N_DISK - 2; i++) {
0111     exactprojdisk(settings_.zmean(settings_.projdisks(iSeed_, i) - 1),
0112                   rinv,
0113                   phi0,
0114                   t,
0115                   z0,
0116                   phiproj[i],
0117                   rproj[i],
0118                   phider[i],
0119                   rder[i]);
0120   }
0121 
0122   for (unsigned int i = 0; i < N_DISK - 2; i++) {
0123     exactproj(settings_.rmean(i), rinv, phi0, t, z0, phiprojLayer[i], zprojLayer[i], phiderLayer[i], zderLayer[i]);
0124   }
0125 }
0126 
0127 void TrackletCalculatorBase::exacttrackletOverlap(double r1,
0128                                                   double z1,
0129                                                   double phi1,
0130                                                   double r2,
0131                                                   double z2,
0132                                                   double phi2,
0133                                                   double,
0134                                                   double& rinv,
0135                                                   double& phi0,
0136                                                   double& t,
0137                                                   double& z0,
0138                                                   double phiprojLayer[N_PSLAYER],
0139                                                   double zprojLayer[N_PSLAYER],
0140                                                   double phiderLayer[N_PSLAYER],
0141                                                   double zderLayer[N_PSLAYER],
0142                                                   double phiproj[N_DISK - 2],
0143                                                   double rproj[N_DISK - 2],
0144                                                   double phider[N_DISK - 2],
0145                                                   double rder[N_DISK - 2]) {
0146   double deltaphi = reco::reduceRange(phi1 - phi2);
0147 
0148   double dist = sqrt(r2 * r2 + r1 * r1 - 2 * r1 * r2 * cos(deltaphi));
0149 
0150   rinv = 2 * sin(deltaphi) / dist;
0151 
0152   if (r1 > r2)
0153     rinv = -rinv;
0154 
0155   double phi1tmp = phi1 - phimin_;
0156 
0157   phi0 = reco::reduceRange(phi1tmp + asin(0.5 * r1 * rinv));
0158 
0159   double rhopsi1 = 2 * asin(0.5 * r1 * rinv) / rinv;
0160   double rhopsi2 = 2 * asin(0.5 * r2 * rinv) / rinv;
0161 
0162   t = (z1 - z2) / (rhopsi1 - rhopsi2);
0163 
0164   z0 = z1 - t * rhopsi1;
0165 
0166   for (int i = 0; i < 4; i++) {
0167     exactprojdisk(settings_.zmean(i + 1), rinv, phi0, t, z0, phiproj[i], rproj[i], phider[i], rder[i]);
0168   }
0169 
0170   for (int i = 0; i < 1; i++) {
0171     exactproj(settings_.rmean(i), rinv, phi0, t, z0, phiprojLayer[i], zprojLayer[i], phiderLayer[i], zderLayer[i]);
0172   }
0173 }
0174 
0175 void TrackletCalculatorBase::exactproj(double rproj,
0176                                        double rinv,
0177                                        double phi0,
0178                                        double t,
0179                                        double z0,
0180                                        double& phiproj,
0181                                        double& zproj,
0182                                        double& phider,
0183                                        double& zder) {
0184   phiproj = phi0 - asin(0.5 * rproj * rinv);
0185   zproj = z0 + (2 * t / rinv) * asin(0.5 * rproj * rinv);
0186 
0187   phider = -0.5 * rinv / sqrt(1 - pow(0.5 * rproj * rinv, 2));
0188   zder = t / sqrt(1 - pow(0.5 * rproj * rinv, 2));
0189 }
0190 
0191 void TrackletCalculatorBase::exactprojdisk(double zproj,
0192                                            double rinv,
0193                                            double phi0,
0194                                            double t,
0195                                            double z0,
0196                                            double& phiproj,
0197                                            double& rproj,
0198                                            double& phider,
0199                                            double& rder) {
0200   if (t < 0)
0201     zproj = -zproj;
0202 
0203   double tmp = rinv * (zproj - z0) / (2.0 * t);
0204   rproj = (2.0 / rinv) * sin(tmp);
0205   phiproj = phi0 - tmp;
0206 
0207   phider = -rinv / (2 * t);
0208   rder = cos(tmp) / t;
0209 }
0210 
0211 void TrackletCalculatorBase::addDiskProj(Tracklet* tracklet, int disk) {
0212   disk = std::abs(disk);
0213 
0214   FPGAWord fpgar = tracklet->proj(N_LAYER + disk - 1).fpgarzproj();
0215 
0216   if (fpgar.value() * settings_.krprojshiftdisk() < settings_.rmindiskvm())
0217     return;
0218   if (fpgar.value() * settings_.krprojshiftdisk() > settings_.rmaxdisk())
0219     return;
0220 
0221   FPGAWord fpgaphi = tracklet->proj(N_LAYER + disk - 1).fpgaphiproj();
0222 
0223   int iphivmRaw = fpgaphi.value() >> (fpgaphi.nbits() - 5);
0224 
0225   int iphi = iphivmRaw / (32 / settings_.nallstubs(disk + N_LAYER - 1));
0226 
0227   addProjectionDisk(disk, iphi, trackletprojdisks_[disk - 1][iphi], tracklet);
0228 }
0229 
0230 bool TrackletCalculatorBase::addLayerProj(Tracklet* tracklet, int layer) {
0231   assert(layer > 0);
0232 
0233   FPGAWord fpgaz = tracklet->proj(layer - 1).fpgarzproj();
0234   FPGAWord fpgaphi = tracklet->proj(layer - 1).fpgaphiproj();
0235 
0236   if (fpgaphi.atExtreme())
0237     edm::LogProblem("Tracklet") << "at extreme! " << fpgaphi.value();
0238 
0239   assert(!fpgaphi.atExtreme());
0240 
0241   if (fpgaz.atExtreme())
0242     return false;
0243 
0244   if (std::abs(fpgaz.value() * settings_.kz()) > settings_.zlength())
0245     return false;
0246 
0247   int iphivmRaw = fpgaphi.value() >> (fpgaphi.nbits() - 5);
0248   int iphi = iphivmRaw / (32 / settings_.nallstubs(layer - 1));
0249 
0250   addProjection(layer, iphi, trackletprojlayers_[layer - 1][iphi], tracklet);
0251 
0252   return true;
0253 }
0254 
0255 void TrackletCalculatorBase::addProjection(int layer,
0256                                            int iphi,
0257                                            TrackletProjectionsMemory* trackletprojs,
0258                                            Tracklet* tracklet) {
0259   if (trackletprojs == nullptr) {
0260     if (settings_.warnNoMem()) {
0261       edm::LogVerbatim("Tracklet") << "No projection memory exists in " << getName() << " for layer = " << layer
0262                                    << " iphi = " << iphi + 1;
0263     }
0264     return;
0265   }
0266   assert(trackletprojs != nullptr);
0267   trackletprojs->addProj(tracklet);
0268 }
0269 
0270 void TrackletCalculatorBase::addProjectionDisk(int disk,
0271                                                int iphi,
0272                                                TrackletProjectionsMemory* trackletprojs,
0273                                                Tracklet* tracklet) {
0274   if (iSeed_ == Seed::L3L4 && abs(disk) == 4)
0275     return;  //L3L4 projections to D3 are not used. Should be in configuration
0276   if (trackletprojs == nullptr) {
0277     if (iSeed_ == Seed::L3L4 && abs(disk) == 3)
0278       return;  //L3L4 projections to D3 are not used.
0279     if (settings_.warnNoMem()) {
0280       edm::LogVerbatim("Tracklet") << "No projection memory exists in " << getName() << " for disk = " << abs(disk)
0281                                    << " iphi = " << iphi + 1;
0282     }
0283     return;
0284   }
0285   assert(trackletprojs != nullptr);
0286   trackletprojs->addProj(tracklet);
0287 }
0288 
0289 bool TrackletCalculatorBase::goodTrackPars(bool goodrinv, bool goodz0) {
0290   bool success = true;
0291   if (!goodrinv) {
0292     if (settings_.debugTracklet()) {
0293       edm::LogVerbatim("Tracklet") << getName() << " TrackletCalculatorBase irinv too large";
0294     }
0295     success = false;
0296   }
0297   if (!goodz0) {
0298     if (settings_.debugTracklet()) {
0299       edm::LogVerbatim("Tracklet") << getName() << " TrackletCalculatorBase z0 cut to large";
0300     }
0301     success = false;
0302   }
0303   return success;
0304 }
0305 
0306 bool TrackletCalculatorBase::inSector(int iphi0, int irinv, double phi0approx, double rinvapprox) {
0307   double phicritapprox = phi0approx - asin(0.5 * settings_.rcrit() * rinvapprox);
0308 
0309   int ifactor = 0.5 * settings_.rcrit() * settings_.krinvpars() / settings_.kphi0pars() * (1 << 8);
0310   int iphicrit = iphi0 - (irinv >> 8) * ifactor;
0311 
0312   int iphicritmincut = settings_.phicritminmc() / globals_->ITC_L1L2()->phi0_final.K();
0313   int iphicritmaxcut = settings_.phicritmaxmc() / globals_->ITC_L1L2()->phi0_final.K();
0314 
0315   bool keepapprox = (phicritapprox > settings_.phicritminmc()) && (phicritapprox < settings_.phicritmaxmc()),
0316        keep = (iphicrit > iphicritmincut) && (iphicrit < iphicritmaxcut);
0317   if (settings_.debugTracklet())
0318     if (keepapprox && !keep)
0319       edm::LogVerbatim("Tracklet") << getName()
0320                                    << " Tracklet kept with exact phicrit cut but not approximate, phicritapprox: "
0321                                    << phicritapprox;
0322   if (settings_.usephicritapprox()) {
0323     return keepapprox;
0324   } else {
0325     return keep;
0326   }
0327 
0328   return true;
0329 }
0330 
0331 bool TrackletCalculatorBase::barrelSeeding(const Stub* innerFPGAStub,
0332                                            const L1TStub* innerStub,
0333                                            const Stub* outerFPGAStub,
0334                                            const L1TStub* outerStub) {
0335   if (settings_.debugTracklet()) {
0336     edm::LogVerbatim("Tracklet") << "TrackletCalculatorBase " << getName()
0337                                  << " trying stub pair in layer (inner outer): " << innerFPGAStub->layer().value()
0338                                  << " " << outerFPGAStub->layer().value();
0339   }
0340 
0341   assert(outerFPGAStub->layerdisk() < N_LAYER);
0342   assert(layerdisk1_ == (unsigned int)innerFPGAStub->layer().value());
0343   assert(layerdisk1_ < N_LAYER && layerdisk2_ < N_LAYER);
0344 
0345   double r1 = innerStub->r();
0346   double z1 = innerStub->z();
0347   double phi1 = innerStub->phi();
0348 
0349   double r2 = outerStub->r();
0350   double z2 = outerStub->z();
0351   double phi2 = outerStub->phi();
0352 
0353   double rinv, phi0, t, z0;
0354 
0355   double phiproj[N_LAYER - 2], zproj[N_LAYER - 2], phider[N_LAYER - 2], zder[N_LAYER - 2];
0356   double phiprojdisk[N_DISK], rprojdisk[N_DISK], phiderdisk[N_DISK], rderdisk[N_DISK];
0357 
0358   exacttracklet(r1,
0359                 z1,
0360                 phi1,
0361                 r2,
0362                 z2,
0363                 phi2,
0364                 outerStub->sigmaz(),
0365                 rinv,
0366                 phi0,
0367                 t,
0368                 z0,
0369                 phiproj,
0370                 zproj,
0371                 phider,
0372                 zder,
0373                 phiprojdisk,
0374                 rprojdisk,
0375                 phiderdisk,
0376                 rderdisk);
0377 
0378   if (settings_.useapprox()) {
0379     phi1 = innerFPGAStub->phiapprox(phimin_, phimax_);
0380     z1 = innerFPGAStub->zapprox();
0381     r1 = innerFPGAStub->rapprox();
0382 
0383     phi2 = outerFPGAStub->phiapprox(phimin_, phimax_);
0384     z2 = outerFPGAStub->zapprox();
0385     r2 = outerFPGAStub->rapprox();
0386   }
0387 
0388   double rinvapprox, phi0approx, tapprox, z0approx;
0389   double phiprojapprox[N_LAYER - 2], zprojapprox[N_LAYER - 2];
0390   double phiprojdiskapprox[N_DISK], rprojdiskapprox[N_DISK];
0391 
0392   IMATH_TrackletCalculator* ITC;
0393   if (iSeed_ == 0)
0394     ITC = globals_->ITC_L1L2();
0395   else if (iSeed_ == 1)
0396     ITC = globals_->ITC_L2L3();
0397   else if (iSeed_ == 2)
0398     ITC = globals_->ITC_L3L4();
0399   else
0400     ITC = globals_->ITC_L5L6();
0401 
0402   ITC->r1.set_fval(r1 - settings_.rmean(layerdisk1_));
0403   ITC->r2.set_fval(r2 - settings_.rmean(layerdisk2_));
0404   ITC->z1.set_fval(z1);
0405   ITC->z2.set_fval(z2);
0406   double sphi1 = angle0to2pi::make0To2pi(phi1 - phimin_);
0407   double sphi2 = angle0to2pi::make0To2pi(phi2 - phimin_);
0408 
0409   ITC->phi1.set_fval(sphi1);
0410   ITC->phi2.set_fval(sphi2);
0411 
0412   ITC->rproj0.set_fval(settings_.rmean(settings_.projlayers(iSeed_, 0) - 1));
0413   ITC->rproj1.set_fval(settings_.rmean(settings_.projlayers(iSeed_, 1) - 1));
0414   ITC->rproj2.set_fval(settings_.rmean(settings_.projlayers(iSeed_, 2) - 1));
0415   ITC->rproj3.set_fval(settings_.rmean(settings_.projlayers(iSeed_, 3) - 1));
0416 
0417   ITC->zproj0.set_fval(t > 0 ? settings_.zmean(0) : -settings_.zmean(0));
0418   ITC->zproj1.set_fval(t > 0 ? settings_.zmean(1) : -settings_.zmean(1));
0419   ITC->zproj2.set_fval(t > 0 ? settings_.zmean(2) : -settings_.zmean(2));
0420   ITC->zproj3.set_fval(t > 0 ? settings_.zmean(3) : -settings_.zmean(3));
0421   ITC->zproj4.set_fval(t > 0 ? settings_.zmean(4) : -settings_.zmean(4));
0422 
0423   ITC->rinv_final.calculate();
0424   ITC->phi0_final.calculate();
0425   ITC->t_final.calculate();
0426   ITC->z0_final.calculate();
0427 
0428   ITC->phiL_0_final.calculate();
0429   ITC->phiL_1_final.calculate();
0430   ITC->phiL_2_final.calculate();
0431   ITC->phiL_3_final.calculate();
0432 
0433   ITC->zL_0_final.calculate();
0434   ITC->zL_1_final.calculate();
0435   ITC->zL_2_final.calculate();
0436   ITC->zL_3_final.calculate();
0437 
0438   ITC->phiD_0_final.calculate();
0439   ITC->phiD_1_final.calculate();
0440   ITC->phiD_2_final.calculate();
0441   ITC->phiD_3_final.calculate();
0442   ITC->phiD_4_final.calculate();
0443 
0444   ITC->rD_0_final.calculate();
0445   ITC->rD_1_final.calculate();
0446   ITC->rD_2_final.calculate();
0447   ITC->rD_3_final.calculate();
0448   ITC->rD_4_final.calculate();
0449 
0450   ITC->der_phiL_final.calculate();
0451   ITC->der_zL_final.calculate();
0452   ITC->der_phiD_final.calculate();
0453   ITC->der_rD_final.calculate();
0454 
0455   //store the approximate results
0456   rinvapprox = ITC->rinv_final.fval();
0457   phi0approx = ITC->phi0_final.fval();
0458   tapprox = ITC->t_final.fval();
0459   z0approx = ITC->z0_final.fval();
0460 
0461   phiprojapprox[0] = ITC->phiL_0_final.fval();
0462   phiprojapprox[1] = ITC->phiL_1_final.fval();
0463   phiprojapprox[2] = ITC->phiL_2_final.fval();
0464   phiprojapprox[3] = ITC->phiL_3_final.fval();
0465 
0466   zprojapprox[0] = ITC->zL_0_final.fval();
0467   zprojapprox[1] = ITC->zL_1_final.fval();
0468   zprojapprox[2] = ITC->zL_2_final.fval();
0469   zprojapprox[3] = ITC->zL_3_final.fval();
0470 
0471   phiprojdiskapprox[0] = ITC->phiD_0_final.fval();
0472   phiprojdiskapprox[1] = ITC->phiD_1_final.fval();
0473   phiprojdiskapprox[2] = ITC->phiD_2_final.fval();
0474   phiprojdiskapprox[3] = ITC->phiD_3_final.fval();
0475   phiprojdiskapprox[4] = ITC->phiD_4_final.fval();
0476 
0477   rprojdiskapprox[0] = ITC->rD_0_final.fval();
0478   rprojdiskapprox[1] = ITC->rD_1_final.fval();
0479   rprojdiskapprox[2] = ITC->rD_2_final.fval();
0480   rprojdiskapprox[3] = ITC->rD_3_final.fval();
0481   rprojdiskapprox[4] = ITC->rD_4_final.fval();
0482 
0483   //now binary
0484 
0485   int irinv, iphi0, it, iz0;
0486   Projection projs[N_LAYER + N_DISK];
0487 
0488   int iphiproj[N_LAYER - 2], izproj[N_LAYER - 2];
0489   int iphiprojdisk[N_DISK], irprojdisk[N_DISK];
0490 
0491   int ir1 = innerFPGAStub->r().value();
0492   int iphi1 = innerFPGAStub->phi().value();
0493   int iz1 = innerFPGAStub->z().value();
0494 
0495   int ir2 = outerFPGAStub->r().value();
0496   int iphi2 = outerFPGAStub->phi().value();
0497   int iz2 = outerFPGAStub->z().value();
0498 
0499   iphi1 <<= (settings_.nphibitsstub(5) - settings_.nphibitsstub(layerdisk1_));
0500   iphi2 <<= (settings_.nphibitsstub(5) - settings_.nphibitsstub(layerdisk2_));
0501   ir1 <<= (8 - settings_.nrbitsstub(layerdisk1_));
0502   ir2 <<= (8 - settings_.nrbitsstub(layerdisk2_));
0503 
0504   iz1 <<= (settings_.nzbitsstub(0) - settings_.nzbitsstub(layerdisk1_));
0505   iz2 <<= (settings_.nzbitsstub(0) - settings_.nzbitsstub(layerdisk2_));
0506 
0507   ITC->r1.set_ival(ir1);
0508   ITC->r2.set_ival(ir2);
0509   ITC->z1.set_ival(iz1);
0510   ITC->z2.set_ival(iz2);
0511   ITC->phi1.set_ival(iphi1);
0512   ITC->phi2.set_ival(iphi2);
0513 
0514   ITC->rinv_final.calculate();
0515   ITC->phi0_final.calculate();
0516   ITC->t_final.calculate();
0517   ITC->z0_final.calculate();
0518 
0519   ITC->phiL_0_final.calculate();
0520   ITC->phiL_1_final.calculate();
0521   ITC->phiL_2_final.calculate();
0522   ITC->phiL_3_final.calculate();
0523 
0524   ITC->zL_0_final.calculate();
0525   ITC->zL_1_final.calculate();
0526   ITC->zL_2_final.calculate();
0527   ITC->zL_3_final.calculate();
0528 
0529   ITC->phiD_0_final.calculate();
0530   ITC->phiD_1_final.calculate();
0531   ITC->phiD_2_final.calculate();
0532   ITC->phiD_3_final.calculate();
0533   ITC->phiD_4_final.calculate();
0534 
0535   ITC->rD_0_final.calculate();
0536   ITC->rD_1_final.calculate();
0537   ITC->rD_2_final.calculate();
0538   ITC->rD_3_final.calculate();
0539   ITC->rD_4_final.calculate();
0540 
0541   ITC->der_phiL_final.calculate();
0542   ITC->der_zL_final.calculate();
0543   ITC->der_phiD_final.calculate();
0544   ITC->der_rD_final.calculate();
0545 
0546   //store the binary results
0547   irinv = ITC->rinv_final.ival();
0548   iphi0 = ITC->phi0_final.ival();
0549   it = ITC->t_final.ival();
0550   iz0 = ITC->z0_final.ival();
0551 
0552   iphiproj[0] = ITC->phiL_0_final.ival();
0553   iphiproj[1] = ITC->phiL_1_final.ival();
0554   iphiproj[2] = ITC->phiL_2_final.ival();
0555   iphiproj[3] = ITC->phiL_3_final.ival();
0556 
0557   izproj[0] = ITC->zL_0_final.ival();
0558   izproj[1] = ITC->zL_1_final.ival();
0559   izproj[2] = ITC->zL_2_final.ival();
0560   izproj[3] = ITC->zL_3_final.ival();
0561 
0562   if (!goodTrackPars(ITC->rinv_final.local_passes(), ITC->z0_final.local_passes())) {
0563     if (settings_.debugTracklet()) {
0564       edm::LogVerbatim("Tracklet") << getName() << " Failed rinv or z0 cut";
0565     }
0566     return false;
0567   }
0568 
0569   if (!inSector(iphi0, irinv, phi0approx, rinvapprox)) {
0570     if (settings_.debugTracklet()) {
0571       edm::LogVerbatim("Tracklet") << getName() << " Failed in sector check";
0572     }
0573     return false;
0574   }
0575 
0576   for (unsigned int i = 0; i < N_LAYER - 2; ++i) {
0577     //reject projection if z is out of range
0578     if (izproj[i] < -(1 << (settings_.nzbitsstub(0) - 1)))
0579       continue;
0580     if (izproj[i] >= (1 << (settings_.nzbitsstub(0) - 1)))
0581       continue;
0582 
0583     //reject projection if phi is out of range
0584     if (iphiproj[i] >= (1 << settings_.nphibitsstub(5)) - 1)
0585       continue;
0586     if (iphiproj[i] <= 0)
0587       continue;
0588 
0589     //Adjust bits for r and z projection depending on layer
0590     if (settings_.projlayers(iSeed_, i) <= 3) {  //TODO clean up logic
0591       iphiproj[i] >>= (settings_.nphibitsstub(5) - settings_.nphibitsstub(settings_.projlayers(iSeed_, i) - 1));
0592     } else {
0593       izproj[i] >>= (settings_.nzbitsstub(0) - settings_.nzbitsstub(5));
0594     }
0595 
0596     projs[settings_.projlayers(iSeed_, i) - 1].init(settings_,
0597                                                     settings_.projlayers(iSeed_, i) - 1,
0598                                                     iphiproj[i],
0599                                                     izproj[i],
0600                                                     ITC->der_phiL_final.ival(),
0601                                                     ITC->der_zL_final.ival(),
0602                                                     phiproj[i],
0603                                                     zproj[i],
0604                                                     phider[i],
0605                                                     zder[i],
0606                                                     phiprojapprox[i],
0607                                                     zprojapprox[i],
0608                                                     ITC->der_phiL_final.fval(),
0609                                                     ITC->der_zL_final.fval(),
0610                                                     !(iSeed_ == 2 || iSeed_ == 3));
0611   }
0612 
0613   iphiprojdisk[0] = ITC->phiD_0_final.ival();
0614   iphiprojdisk[1] = ITC->phiD_1_final.ival();
0615   iphiprojdisk[2] = ITC->phiD_2_final.ival();
0616   iphiprojdisk[3] = ITC->phiD_3_final.ival();
0617   iphiprojdisk[4] = ITC->phiD_4_final.ival();
0618 
0619   irprojdisk[0] = ITC->rD_0_final.ival();
0620   irprojdisk[1] = ITC->rD_1_final.ival();
0621   irprojdisk[2] = ITC->rD_2_final.ival();
0622   irprojdisk[3] = ITC->rD_3_final.ival();
0623   irprojdisk[4] = ITC->rD_4_final.ival();
0624 
0625   if (std::abs(it * ITC->t_final.K()) > 1.0) {
0626     for (unsigned int i = 0; i < N_DISK; ++i) {
0627       if (iphiprojdisk[i] <= 0)
0628         continue;
0629       if (iphiprojdisk[i] >= (1 << settings_.nphibitsstub(0)) - 1)
0630         continue;
0631 
0632       if (irprojdisk[i] < settings_.rmindisk() / ITC->rD_0_final.K() ||
0633           irprojdisk[i] >= settings_.rmaxdisk() / ITC->rD_0_final.K())
0634         continue;
0635 
0636       projs[i + N_LAYER].init(settings_,
0637                               i + N_LAYER,
0638                               iphiprojdisk[i],
0639                               irprojdisk[i],
0640                               ITC->der_phiD_final.ival(),
0641                               ITC->der_rD_final.ival(),
0642                               phiprojdisk[i],
0643                               rprojdisk[i],
0644                               phiderdisk[i],
0645                               rderdisk[i],
0646                               phiprojdiskapprox[i],
0647                               rprojdiskapprox[i],
0648                               ITC->der_phiD_final.fval(),
0649                               ITC->der_rD_final.fval(),
0650                               !(iSeed_ == 2 || iSeed_ == 3));
0651     }
0652   }
0653 
0654   if (settings_.writeMonitorData("TPars")) {
0655     globals_->ofstream("trackletpars.txt")
0656         << "Trackpars " << layerdisk1_ + 1 << "   " << rinv << " " << rinvapprox << " " << ITC->rinv_final.fval()
0657         << "   " << phi0 << " " << phi0approx << " " << ITC->phi0_final.fval() << "   " << t << " " << tapprox << " "
0658         << ITC->t_final.fval() << "   " << z0 << " " << z0approx << " " << ITC->z0_final.fval() << endl;
0659   }
0660 
0661   Tracklet* tracklet = new Tracklet(settings_,
0662                                     iSeed_,
0663                                     innerFPGAStub,
0664                                     nullptr,
0665                                     outerFPGAStub,
0666                                     rinv,
0667                                     phi0,
0668                                     0.0,
0669                                     z0,
0670                                     t,
0671                                     rinvapprox,
0672                                     phi0approx,
0673                                     0.0,
0674                                     z0approx,
0675                                     tapprox,
0676                                     irinv,
0677                                     iphi0,
0678                                     0,
0679                                     iz0,
0680                                     it,
0681                                     projs,
0682                                     false);
0683 
0684   if (settings_.debugTracklet()) {
0685     edm::LogVerbatim("Tracklet") << "TrackletCalculator " << getName() << " Found tracklet for seed = " << iSeed_ << " "
0686                                  << iSector_ << " phi0 = " << phi0;
0687   }
0688 
0689   tracklet->setTrackletIndex(trackletpars_->nTracklets());
0690   tracklet->setTCIndex(TCIndex_);
0691 
0692   if (settings_.writeMonitorData("Seeds")) {
0693     ofstream fout("seeds.txt", ofstream::app);
0694     fout << __FILE__ << ":" << __LINE__ << " " << name_ << "_" << iSector_ << " " << tracklet->getISeed() << endl;
0695     fout.close();
0696   }
0697   trackletpars_->addTracklet(tracklet);
0698 
0699   if (settings_.bookHistos()) {
0700     HistBase* hists = globals_->histograms();
0701     int tp = tracklet->tpseed();
0702     hists->fillTrackletParams(settings_,
0703                               globals_,
0704                               iSeed_,
0705                               iSector_,
0706                               rinvapprox,
0707                               irinv * ITC->rinv_final.K(),
0708                               phi0approx,
0709                               iphi0 * ITC->phi0_final.K(),
0710                               asinh(tapprox),
0711                               asinh(it * ITC->t_final.K()),
0712                               z0approx,
0713                               iz0 * ITC->z0_final.K(),
0714                               tp);
0715   }
0716 
0717   bool addL3 = false;
0718   bool addL4 = false;
0719   bool addL5 = false;
0720   bool addL6 = false;
0721   for (unsigned int j = 0; j < N_LAYER - 2; j++) {
0722     int lproj = settings_.projlayers(iSeed_, j);
0723     bool added = false;
0724     if (tracklet->validProj(lproj - 1)) {
0725       added = addLayerProj(tracklet, lproj);
0726       if (added && lproj == 3)
0727         addL3 = true;
0728       if (added && lproj == 4)
0729         addL4 = true;
0730       if (added && lproj == 5)
0731         addL5 = true;
0732       if (added && lproj == 6)
0733         addL6 = true;
0734     }
0735   }
0736 
0737   for (unsigned int j = 0; j < N_DISK - 1; j++) {  //no projections to 5th disk!!
0738     int disk = j + 1;
0739     if (disk == 4 && addL3)
0740       continue;
0741     if (disk == 3 && addL4)
0742       continue;
0743     if (disk == 2 && addL5)
0744       continue;
0745     if (disk == 1 && addL6)
0746       continue;
0747     if (it < 0)
0748       disk = -disk;
0749     if (tracklet->validProj(N_LAYER + abs(disk) - 1)) {
0750       addDiskProj(tracklet, disk);
0751     }
0752   }
0753 
0754   return true;
0755 }
0756 
0757 bool TrackletCalculatorBase::diskSeeding(const Stub* innerFPGAStub,
0758                                          const L1TStub* innerStub,
0759                                          const Stub* outerFPGAStub,
0760                                          const L1TStub* outerStub) {
0761   if (settings_.debugTracklet()) {
0762     edm::LogVerbatim("Tracklet") << "TrackletCalculator::execute calculate disk seeds";
0763   }
0764 
0765   int sign = 1;
0766   if (innerFPGAStub->disk().value() < 0)
0767     sign = -1;
0768 
0769   int disk = innerFPGAStub->disk().value();
0770   assert(abs(disk) == 1 || abs(disk) == 3);
0771 
0772   assert(innerStub->isPSmodule());
0773   assert(outerStub->isPSmodule());
0774 
0775   double r1 = innerStub->r();
0776   double z1 = innerStub->z();
0777   double phi1 = innerStub->phi();
0778 
0779   double r2 = outerStub->r();
0780   double z2 = outerStub->z();
0781   double phi2 = outerStub->phi();
0782 
0783   if (r2 < r1 + 2.0) {
0784     return false;  //Protection... Should be handled cleaner to avoid problem with floating point calculation
0785   }
0786 
0787   double rinv, phi0, t, z0;
0788 
0789   double phiproj[N_PSLAYER], zproj[N_PSLAYER], phider[N_PSLAYER], zder[N_PSLAYER];
0790   double phiprojdisk[N_DISK - 2], rprojdisk[N_DISK - 2], phiderdisk[N_DISK - 2], rderdisk[N_DISK - 2];
0791 
0792   exacttrackletdisk(r1,
0793                     z1,
0794                     phi1,
0795                     r2,
0796                     z2,
0797                     phi2,
0798                     outerStub->sigmaz(),
0799                     rinv,
0800                     phi0,
0801                     t,
0802                     z0,
0803                     phiproj,
0804                     zproj,
0805                     phider,
0806                     zder,
0807                     phiprojdisk,
0808                     rprojdisk,
0809                     phiderdisk,
0810                     rderdisk);
0811 
0812   //Truncates floating point positions to integer representation precision
0813   if (settings_.useapprox()) {
0814     phi1 = innerFPGAStub->phiapprox(phimin_, phimax_);
0815     z1 = innerFPGAStub->zapprox();
0816     r1 = innerFPGAStub->rapprox();
0817 
0818     phi2 = outerFPGAStub->phiapprox(phimin_, phimax_);
0819     z2 = outerFPGAStub->zapprox();
0820     r2 = outerFPGAStub->rapprox();
0821   }
0822 
0823   double rinvapprox, phi0approx, tapprox, z0approx;
0824   double phiprojapprox[N_PSLAYER], zprojapprox[N_PSLAYER];
0825   double phiprojdiskapprox[N_DISK - 2], rprojdiskapprox[N_DISK - 2];
0826 
0827   IMATH_TrackletCalculatorDisk* ITC;
0828   if (disk == 1)
0829     ITC = globals_->ITC_F1F2();
0830   else if (disk == 3)
0831     ITC = globals_->ITC_F3F4();
0832   else if (disk == -1)
0833     ITC = globals_->ITC_B1B2();
0834   else
0835     ITC = globals_->ITC_B3B4();
0836 
0837   ITC->r1.set_fval(r1);
0838   ITC->r2.set_fval(r2);
0839   int signt = t > 0 ? 1 : -1;
0840   ITC->z1.set_fval(z1 - signt * settings_.zmean(layerdisk1_ - N_LAYER));
0841   ITC->z2.set_fval(z2 - signt * settings_.zmean(layerdisk2_ - N_LAYER));
0842   double sphi1 = angle0to2pi::make0To2pi(phi1 - phimin_);
0843   double sphi2 = angle0to2pi::make0To2pi(phi2 - phimin_);
0844   ITC->phi1.set_fval(sphi1);
0845   ITC->phi2.set_fval(sphi2);
0846 
0847   ITC->rproj0.set_fval(settings_.rmean(0));
0848   ITC->rproj1.set_fval(settings_.rmean(1));
0849   ITC->rproj2.set_fval(settings_.rmean(2));
0850 
0851   ITC->zproj0.set_fval(signt * settings_.zmean(settings_.projdisks(iSeed_, 0) - 1));
0852   ITC->zproj1.set_fval(signt * settings_.zmean(settings_.projdisks(iSeed_, 1) - 1));
0853   ITC->zproj2.set_fval(signt * settings_.zmean(settings_.projdisks(iSeed_, 2) - 1));
0854 
0855   ITC->rinv_final.calculate();
0856   ITC->phi0_final.calculate();
0857   ITC->t_final.calculate();
0858   ITC->z0_final.calculate();
0859 
0860   ITC->phiL_0_final.calculate();
0861   ITC->phiL_1_final.calculate();
0862   ITC->phiL_2_final.calculate();
0863 
0864   ITC->zL_0_final.calculate();
0865   ITC->zL_1_final.calculate();
0866   ITC->zL_2_final.calculate();
0867 
0868   ITC->phiD_0_final.calculate();
0869   ITC->phiD_1_final.calculate();
0870   ITC->phiD_2_final.calculate();
0871 
0872   ITC->rD_0_final.calculate();
0873   ITC->rD_1_final.calculate();
0874   ITC->rD_2_final.calculate();
0875 
0876   ITC->der_phiL_final.calculate();
0877   ITC->der_zL_final.calculate();
0878   ITC->der_phiD_final.calculate();
0879   ITC->der_rD_final.calculate();
0880 
0881   //store the approximate results
0882   rinvapprox = ITC->rinv_final.fval();
0883   phi0approx = ITC->phi0_final.fval();
0884   tapprox = ITC->t_final.fval();
0885   z0approx = ITC->z0_final.fval();
0886 
0887   phiprojapprox[0] = ITC->phiL_0_final.fval();
0888   phiprojapprox[1] = ITC->phiL_1_final.fval();
0889   phiprojapprox[2] = ITC->phiL_2_final.fval();
0890 
0891   zprojapprox[0] = ITC->zL_0_final.fval();
0892   zprojapprox[1] = ITC->zL_1_final.fval();
0893   zprojapprox[2] = ITC->zL_2_final.fval();
0894 
0895   phiprojdiskapprox[0] = ITC->phiD_0_final.fval();
0896   phiprojdiskapprox[1] = ITC->phiD_1_final.fval();
0897   phiprojdiskapprox[2] = ITC->phiD_2_final.fval();
0898 
0899   rprojdiskapprox[0] = ITC->rD_0_final.fval();
0900   rprojdiskapprox[1] = ITC->rD_1_final.fval();
0901   rprojdiskapprox[2] = ITC->rD_2_final.fval();
0902 
0903   //now binary
0904 
0905   int irinv, iphi0, it, iz0;
0906   int iphiproj[N_PSLAYER], izproj[N_PSLAYER];
0907 
0908   int iphiprojdisk[N_DISK - 2], irprojdisk[N_DISK - 2];
0909 
0910   int ir1 = innerFPGAStub->r().value();
0911   int iphi1 = innerFPGAStub->phi().value();
0912   int iz1 = innerFPGAStub->z().value();
0913 
0914   int ir2 = outerFPGAStub->r().value();
0915   int iphi2 = outerFPGAStub->phi().value();
0916   int iz2 = outerFPGAStub->z().value();
0917 
0918   //To get same precission as for layers.
0919   iphi1 <<= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
0920   iphi2 <<= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
0921 
0922   ITC->r1.set_ival(ir1);
0923   ITC->r2.set_ival(ir2);
0924   ITC->z1.set_ival(iz1);
0925   ITC->z2.set_ival(iz2);
0926   ITC->phi1.set_ival(iphi1);
0927   ITC->phi2.set_ival(iphi2);
0928 
0929   ITC->rinv_final.calculate();
0930   ITC->phi0_final.calculate();
0931   ITC->t_final.calculate();
0932   ITC->z0_final.calculate();
0933 
0934   ITC->phiL_0_final.calculate();
0935   ITC->phiL_1_final.calculate();
0936   ITC->phiL_2_final.calculate();
0937 
0938   ITC->zL_0_final.calculate();
0939   ITC->zL_1_final.calculate();
0940   ITC->zL_2_final.calculate();
0941 
0942   ITC->phiD_0_final.calculate();
0943   ITC->phiD_1_final.calculate();
0944   ITC->phiD_2_final.calculate();
0945 
0946   ITC->rD_0_final.calculate();
0947   ITC->rD_1_final.calculate();
0948   ITC->rD_2_final.calculate();
0949 
0950   ITC->der_phiL_final.calculate();
0951   ITC->der_zL_final.calculate();
0952   ITC->der_phiD_final.calculate();
0953   ITC->der_rD_final.calculate();
0954 
0955   //store the binary results
0956   irinv = ITC->rinv_final.ival();
0957   iphi0 = ITC->phi0_final.ival();
0958   it = ITC->t_final.ival();
0959   iz0 = ITC->z0_final.ival();
0960 
0961   iphiproj[0] = ITC->phiL_0_final.ival();
0962   iphiproj[1] = ITC->phiL_1_final.ival();
0963   iphiproj[2] = ITC->phiL_2_final.ival();
0964 
0965   izproj[0] = ITC->zL_0_final.ival();
0966   izproj[1] = ITC->zL_1_final.ival();
0967   izproj[2] = ITC->zL_2_final.ival();
0968 
0969   if (!goodTrackPars(ITC->rinv_final.local_passes(), ITC->z0_final.local_passes()))
0970     return false;
0971 
0972   if (!inSector(iphi0, irinv, phi0approx, rinvapprox))
0973     return false;
0974 
0975   Projection projs[N_LAYER + N_DISK];
0976 
0977   for (unsigned int i = 0; i < N_DISK - 2; ++i) {
0978     //Check is outside z range
0979     if (izproj[i] < -(1 << (settings_.nzbitsstub(0) - 1)))
0980       continue;
0981     if (izproj[i] >= (1 << (settings_.nzbitsstub(0) - 1)))
0982       continue;
0983 
0984     //Check if outside phi range
0985     if (iphiproj[i] >= (1 << settings_.nphibitsstub(5)) - 1)
0986       continue;
0987     if (iphiproj[i] <= 0)
0988       continue;
0989 
0990     //shift bits - allways in PS modules for disk seeding
0991     iphiproj[i] >>= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
0992 
0993     projs[i].init(settings_,
0994                   i,
0995                   iphiproj[i],
0996                   izproj[i],
0997                   ITC->der_phiL_final.ival(),
0998                   ITC->der_zL_final.ival(),
0999                   phiproj[i],
1000                   zproj[i],
1001                   phider[i],
1002                   zder[i],
1003                   phiprojapprox[i],
1004                   zprojapprox[i],
1005                   ITC->der_phiL_final.fval(),
1006                   ITC->der_zL_final.fval(),
1007                   true);
1008   }
1009 
1010   iphiprojdisk[0] = ITC->phiD_0_final.ival();
1011   iphiprojdisk[1] = ITC->phiD_1_final.ival();
1012   iphiprojdisk[2] = ITC->phiD_2_final.ival();
1013 
1014   irprojdisk[0] = ITC->rD_0_final.ival();
1015   irprojdisk[1] = ITC->rD_1_final.ival();
1016   irprojdisk[2] = ITC->rD_2_final.ival();
1017 
1018   for (unsigned int i = 0; i < N_DISK - 2; ++i) {
1019     //check that phi projection in range
1020     if (iphiprojdisk[i] <= 0)
1021       continue;
1022     if (iphiprojdisk[i] >= (1 << settings_.nphibitsstub(0)) - 1)
1023       continue;
1024 
1025     //check that r projection in range
1026     if (irprojdisk[i] <= 0 || irprojdisk[i] >= settings_.rmaxdisk() / ITC->rD_0_final.K())
1027       continue;
1028 
1029     projs[settings_.projdisks(iSeed_, i) + N_LAYER - 1].init(settings_,
1030                                                              settings_.projdisks(iSeed_, i) + N_LAYER - 1,
1031                                                              iphiprojdisk[i],
1032                                                              irprojdisk[i],
1033                                                              ITC->der_phiD_final.ival(),
1034                                                              ITC->der_rD_final.ival(),
1035                                                              phiprojdisk[i],
1036                                                              rprojdisk[i],
1037                                                              phiderdisk[i],
1038                                                              rderdisk[i],
1039                                                              phiprojdiskapprox[i],
1040                                                              rprojdiskapprox[i],
1041                                                              ITC->der_phiD_final.fval(),
1042                                                              ITC->der_rD_final.fval(),
1043                                                              true);
1044   }
1045 
1046   if (settings_.writeMonitorData("TPars")) {
1047     globals_->ofstream("trackletparsdisk.txt")
1048         << "Trackpars         " << layerdisk1_ - 5 << "   " << rinv << " " << rinvapprox << " "
1049         << ITC->rinv_final.fval() << "   " << phi0 << " " << phi0approx << " " << ITC->phi0_final.fval() << "   " << t
1050         << " " << tapprox << " " << ITC->t_final.fval() << "   " << z0 << " " << z0approx << " " << ITC->z0_final.fval()
1051         << endl;
1052   }
1053 
1054   Tracklet* tracklet = new Tracklet(settings_,
1055                                     iSeed_,
1056                                     innerFPGAStub,
1057                                     nullptr,
1058                                     outerFPGAStub,
1059                                     rinv,
1060                                     phi0,
1061                                     0.0,
1062                                     z0,
1063                                     t,
1064                                     rinvapprox,
1065                                     phi0approx,
1066                                     0.0,
1067                                     z0approx,
1068                                     tapprox,
1069                                     irinv,
1070                                     iphi0,
1071                                     0,
1072                                     iz0,
1073                                     it,
1074                                     projs,
1075                                     true);
1076 
1077   if (settings_.debugTracklet()) {
1078     edm::LogVerbatim("Tracklet") << "Found tracklet for disk seed = " << iSeed_ << " " << tracklet << " " << iSector_;
1079   }
1080 
1081   tracklet->setTrackletIndex(trackletpars_->nTracklets());
1082   tracklet->setTCIndex(TCIndex_);
1083 
1084   if (settings_.writeMonitorData("Seeds")) {
1085     ofstream fout("seeds.txt", ofstream::app);
1086     fout << __FILE__ << ":" << __LINE__ << " " << name_ << "_" << iSector_ << " " << tracklet->getISeed() << endl;
1087     fout.close();
1088   }
1089   trackletpars_->addTracklet(tracklet);
1090 
1091   if (tracklet->validProj(0)) {
1092     addLayerProj(tracklet, 1);
1093   }
1094 
1095   if (tracklet->validProj(1)) {
1096     addLayerProj(tracklet, 2);
1097   }
1098 
1099   for (unsigned int j = 0; j < N_DISK - 2; j++) {
1100     if (tracklet->validProj(N_LAYER + settings_.projdisks(iSeed_, j) - 1)) {
1101       addDiskProj(tracklet, sign * settings_.projdisks(iSeed_, j));
1102     }
1103   }
1104 
1105   return true;
1106 }
1107 
1108 bool TrackletCalculatorBase::overlapSeeding(const Stub* innerFPGAStub,
1109                                             const L1TStub* innerStub,
1110                                             const Stub* outerFPGAStub,
1111                                             const L1TStub* outerStub) {
1112   //Deal with overlap stubs here
1113   assert(outerFPGAStub->layerdisk() < N_LAYER);
1114 
1115   assert(innerFPGAStub->layerdisk() >= N_LAYER);
1116 
1117   int disk = innerFPGAStub->disk().value();
1118 
1119   if (settings_.debugTracklet()) {
1120     edm::LogVerbatim("Tracklet") << "trying to make overlap tracklet for seed = " << iSeed_ << " " << getName();
1121   }
1122 
1123   double r1 = innerStub->r();
1124   double z1 = innerStub->z();
1125   double phi1 = innerStub->phi();
1126 
1127   double r2 = outerStub->r();
1128   double z2 = outerStub->z();
1129   double phi2 = outerStub->phi();
1130 
1131   //Protection for wrong radii. Could be handled cleaner to avoid problem with floating point calculation and with overflows in the integer calculation.
1132   if (r1 < r2 + 1.5) {
1133     return false;
1134   }
1135 
1136   double rinv, phi0, t, z0;
1137 
1138   double phiproj[N_PSLAYER], zproj[N_PSLAYER], phider[N_PSLAYER], zder[N_PSLAYER];
1139   double phiprojdisk[N_DISK - 1], rprojdisk[N_DISK - 1], phiderdisk[N_DISK - 1], rderdisk[N_DISK - 1];
1140 
1141   exacttrackletOverlap(r1,
1142                        z1,
1143                        phi1,
1144                        r2,
1145                        z2,
1146                        phi2,
1147                        outerStub->sigmaz(),
1148                        rinv,
1149                        phi0,
1150                        t,
1151                        z0,
1152                        phiproj,
1153                        zproj,
1154                        phider,
1155                        zder,
1156                        phiprojdisk,
1157                        rprojdisk,
1158                        phiderdisk,
1159                        rderdisk);
1160 
1161   //Truncates floating point positions to integer representation precision
1162   if (settings_.useapprox()) {
1163     phi1 = innerFPGAStub->phiapprox(phimin_, phimax_);
1164     z1 = innerFPGAStub->zapprox();
1165     r1 = innerFPGAStub->rapprox();
1166 
1167     phi2 = outerFPGAStub->phiapprox(phimin_, phimax_);
1168     z2 = outerFPGAStub->zapprox();
1169     r2 = outerFPGAStub->rapprox();
1170   }
1171 
1172   double rinvapprox, phi0approx, tapprox, z0approx;
1173   double phiprojapprox[N_PSLAYER], zprojapprox[N_PSLAYER];
1174   double phiprojdiskapprox[N_DISK - 1], rprojdiskapprox[N_DISK - 1];
1175 
1176   IMATH_TrackletCalculatorOverlap* ITC;
1177   int ll = outerFPGAStub->layer().value() + 1;
1178   if (ll == 1 && disk == 1)
1179     ITC = globals_->ITC_L1F1();
1180   else if (ll == 2 && disk == 1)
1181     ITC = globals_->ITC_L2F1();
1182   else if (ll == 1 && disk == -1)
1183     ITC = globals_->ITC_L1B1();
1184   else if (ll == 2 && disk == -1)
1185     ITC = globals_->ITC_L2B1();
1186   else
1187     throw cms::Exception("LogicError") << __FILE__ << " " << __LINE__ << " Invalid seeding!";
1188 
1189   ITC->r1.set_fval(r2 - settings_.rmean(ll - 1));
1190   ITC->r2.set_fval(r1);
1191   int signt = t > 0 ? 1 : -1;
1192   ITC->z1.set_fval(z2);
1193   ITC->z2.set_fval(z1 - signt * settings_.zmean(layerdisk2_ - N_LAYER));
1194   double sphi1 = angle0to2pi::make0To2pi(phi1 - phimin_);
1195   double sphi2 = angle0to2pi::make0To2pi(phi2 - phimin_);
1196   ITC->phi1.set_fval(sphi2);
1197   ITC->phi2.set_fval(sphi1);
1198 
1199   ITC->rproj0.set_fval(settings_.rmean(0));
1200   ITC->rproj1.set_fval(settings_.rmean(1));
1201   ITC->rproj2.set_fval(settings_.rmean(2));
1202 
1203   ITC->zproj0.set_fval(signt * settings_.zmean(1));
1204   ITC->zproj1.set_fval(signt * settings_.zmean(2));
1205   ITC->zproj2.set_fval(signt * settings_.zmean(3));
1206   ITC->zproj3.set_fval(signt * settings_.zmean(4));
1207 
1208   ITC->rinv_final.calculate();
1209   ITC->phi0_final.calculate();
1210   ITC->t_final.calculate();
1211   ITC->z0_final.calculate();
1212 
1213   ITC->phiL_0_final.calculate();
1214   ITC->phiL_1_final.calculate();
1215   ITC->phiL_2_final.calculate();
1216 
1217   ITC->zL_0_final.calculate();
1218   ITC->zL_1_final.calculate();
1219   ITC->zL_2_final.calculate();
1220 
1221   ITC->phiD_0_final.calculate();
1222   ITC->phiD_1_final.calculate();
1223   ITC->phiD_2_final.calculate();
1224   ITC->phiD_3_final.calculate();
1225 
1226   ITC->rD_0_final.calculate();
1227   ITC->rD_1_final.calculate();
1228   ITC->rD_2_final.calculate();
1229   ITC->rD_3_final.calculate();
1230 
1231   ITC->der_phiL_final.calculate();
1232   ITC->der_zL_final.calculate();
1233   ITC->der_phiD_final.calculate();
1234   ITC->der_rD_final.calculate();
1235 
1236   //store the approximate results
1237   rinvapprox = ITC->rinv_final.fval();
1238   phi0approx = ITC->phi0_final.fval();
1239   tapprox = ITC->t_final.fval();
1240   z0approx = ITC->z0_final.fval();
1241 
1242   phiprojapprox[0] = ITC->phiL_0_final.fval();
1243   phiprojapprox[1] = ITC->phiL_1_final.fval();
1244   phiprojapprox[2] = ITC->phiL_2_final.fval();
1245 
1246   zprojapprox[0] = ITC->zL_0_final.fval();
1247   zprojapprox[1] = ITC->zL_1_final.fval();
1248   zprojapprox[2] = ITC->zL_2_final.fval();
1249 
1250   phiprojdiskapprox[0] = ITC->phiD_0_final.fval();
1251   phiprojdiskapprox[1] = ITC->phiD_1_final.fval();
1252   phiprojdiskapprox[2] = ITC->phiD_2_final.fval();
1253   phiprojdiskapprox[3] = ITC->phiD_3_final.fval();
1254 
1255   rprojdiskapprox[0] = ITC->rD_0_final.fval();
1256   rprojdiskapprox[1] = ITC->rD_1_final.fval();
1257   rprojdiskapprox[2] = ITC->rD_2_final.fval();
1258   rprojdiskapprox[3] = ITC->rD_3_final.fval();
1259 
1260   //now binary
1261 
1262   int irinv, iphi0, it, iz0;
1263   int iphiproj[N_LAYER - 2], izproj[N_LAYER - 2];
1264   int iphiprojdisk[N_DISK], irprojdisk[N_DISK];
1265 
1266   int ir2 = innerFPGAStub->r().value();
1267   int iphi2 = innerFPGAStub->phi().value();
1268   int iz2 = innerFPGAStub->z().value();
1269 
1270   int ir1 = outerFPGAStub->r().value();
1271   int iphi1 = outerFPGAStub->phi().value();
1272   int iz1 = outerFPGAStub->z().value();
1273 
1274   //To get global precission
1275   ir1 = l1t::bitShift(ir1, (8 - settings_.nrbitsstub(ll - 1)));
1276   iphi1 <<= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
1277   iphi2 <<= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
1278 
1279   ITC->r1.set_ival(ir1);
1280   ITC->r2.set_ival(ir2);
1281   ITC->z1.set_ival(iz1);
1282   ITC->z2.set_ival(iz2);
1283   ITC->phi1.set_ival(iphi1);
1284   ITC->phi2.set_ival(iphi2);
1285 
1286   ITC->rinv_final.calculate();
1287   ITC->phi0_final.calculate();
1288   ITC->t_final.calculate();
1289   ITC->z0_final.calculate();
1290 
1291   ITC->phiL_0_final.calculate();
1292   ITC->phiL_1_final.calculate();
1293   ITC->phiL_2_final.calculate();
1294 
1295   ITC->zL_0_final.calculate();
1296   ITC->zL_1_final.calculate();
1297   ITC->zL_2_final.calculate();
1298 
1299   ITC->phiD_0_final.calculate();
1300   ITC->phiD_1_final.calculate();
1301   ITC->phiD_2_final.calculate();
1302   ITC->phiD_3_final.calculate();
1303 
1304   ITC->rD_0_final.calculate();
1305   ITC->rD_1_final.calculate();
1306   ITC->rD_2_final.calculate();
1307   ITC->rD_3_final.calculate();
1308 
1309   ITC->der_phiL_final.calculate();
1310   ITC->der_zL_final.calculate();
1311   ITC->der_phiD_final.calculate();
1312   ITC->der_rD_final.calculate();
1313 
1314   //store the binary results
1315   irinv = ITC->rinv_final.ival();
1316   iphi0 = ITC->phi0_final.ival();
1317   it = ITC->t_final.ival();
1318   iz0 = ITC->z0_final.ival();
1319 
1320   iphiproj[0] = ITC->phiL_0_final.ival();
1321   iphiproj[1] = ITC->phiL_1_final.ival();
1322   iphiproj[2] = ITC->phiL_2_final.ival();
1323 
1324   izproj[0] = ITC->zL_0_final.ival();
1325   izproj[1] = ITC->zL_1_final.ival();
1326   izproj[2] = ITC->zL_2_final.ival();
1327 
1328   iphiprojdisk[0] = ITC->phiD_0_final.ival();
1329   iphiprojdisk[1] = ITC->phiD_1_final.ival();
1330   iphiprojdisk[2] = ITC->phiD_2_final.ival();
1331   iphiprojdisk[3] = ITC->phiD_3_final.ival();
1332 
1333   irprojdisk[0] = ITC->rD_0_final.ival();
1334   irprojdisk[1] = ITC->rD_1_final.ival();
1335   irprojdisk[2] = ITC->rD_2_final.ival();
1336   irprojdisk[3] = ITC->rD_3_final.ival();
1337 
1338   if (!goodTrackPars(ITC->rinv_final.local_passes(), ITC->z0_final.local_passes()))
1339     return false;
1340 
1341   if (!inSector(iphi0, irinv, phi0approx, rinvapprox))
1342     return false;
1343 
1344   Projection projs[N_LAYER + N_DISK];
1345 
1346   for (unsigned int i = 0; i < N_DISK - 2; ++i) {
1347     //check that zproj is in range
1348     if (izproj[i] < -(1 << (settings_.nzbitsstub(0) - 1)))
1349       continue;
1350     if (izproj[i] >= (1 << (settings_.nzbitsstub(0) - 1)))
1351       continue;
1352 
1353     //check that phiproj is in range
1354     if (iphiproj[i] >= (1 << settings_.nphibitsstub(5)) - 1)
1355       continue;
1356     if (iphiproj[i] <= 0)
1357       continue;
1358 
1359     //adjust bits for PS modules (no 2S modules in overlap seeds)
1360     iphiproj[i] >>= (settings_.nphibitsstub(5) - settings_.nphibitsstub(0));
1361 
1362     projs[i].init(settings_,
1363                   i,
1364                   iphiproj[i],
1365                   izproj[i],
1366                   ITC->der_phiL_final.ival(),
1367                   ITC->der_zL_final.ival(),
1368                   phiproj[i],
1369                   zproj[i],
1370                   phider[i],
1371                   zder[i],
1372                   phiprojapprox[i],
1373                   zprojapprox[i],
1374                   ITC->der_phiL_final.fval(),
1375                   ITC->der_zL_final.fval(),
1376                   true);
1377   }
1378 
1379   for (int i = 0; i < 4; ++i) {
1380     //check that phi projection in range
1381     if (iphiprojdisk[i] <= 0)
1382       continue;
1383     if (iphiprojdisk[i] >= (1 << settings_.nphibitsstub(0)) - 1)
1384       continue;
1385 
1386     //check that r projection in range
1387     if (irprojdisk[i] <= 0 || irprojdisk[i] >= settings_.rmaxdisk() / ITC->rD_0_final.K())
1388       continue;
1389 
1390     projs[N_LAYER + i + 1].init(settings_,
1391                                 N_LAYER + i + 1,
1392                                 iphiprojdisk[i],
1393                                 irprojdisk[i],
1394                                 ITC->der_phiD_final.ival(),
1395                                 ITC->der_rD_final.ival(),
1396                                 phiprojdisk[i],
1397                                 rprojdisk[i],
1398                                 phiderdisk[i],
1399                                 rderdisk[i],
1400                                 phiprojdiskapprox[i],
1401                                 rprojdiskapprox[i],
1402                                 ITC->der_phiD_final.fval(),
1403                                 ITC->der_rD_final.fval(),
1404                                 true);
1405   }
1406 
1407   if (settings_.writeMonitorData("TPars")) {
1408     globals_->ofstream("trackletparsoverlap.txt")
1409         << "Trackpars " << layerdisk1_ - 5 << "   " << rinv << " " << irinv << " " << ITC->rinv_final.fval() << "   "
1410         << phi0 << " " << iphi0 << " " << ITC->phi0_final.fval() << "   " << t << " " << it << " "
1411         << ITC->t_final.fval() << "   " << z0 << " " << iz0 << " " << ITC->z0_final.fval() << endl;
1412   }
1413 
1414   Tracklet* tracklet = new Tracklet(settings_,
1415                                     iSeed_,
1416                                     innerFPGAStub,
1417                                     nullptr,
1418                                     outerFPGAStub,
1419                                     rinv,
1420                                     phi0,
1421                                     0.0,
1422                                     z0,
1423                                     t,
1424                                     rinvapprox,
1425                                     phi0approx,
1426                                     0.0,
1427                                     z0approx,
1428                                     tapprox,
1429                                     irinv,
1430                                     iphi0,
1431                                     0,
1432                                     iz0,
1433                                     it,
1434                                     projs,
1435                                     false,
1436                                     true);
1437 
1438   if (settings_.debugTracklet()) {
1439     edm::LogVerbatim("Tracklet") << "Found tracklet in overlap seed = " << iSeed_ << " " << tracklet << " " << iSector_;
1440   }
1441 
1442   tracklet->setTrackletIndex(trackletpars_->nTracklets());
1443   tracklet->setTCIndex(TCIndex_);
1444 
1445   if (settings_.writeMonitorData("Seeds")) {
1446     ofstream fout("seeds.txt", ofstream::app);
1447     fout << __FILE__ << ":" << __LINE__ << " " << name_ << "_" << iSector_ << " " << tracklet->getISeed() << endl;
1448     fout.close();
1449   }
1450   trackletpars_->addTracklet(tracklet);
1451 
1452   int layer = outerFPGAStub->layer().value() + 1;
1453 
1454   if (layer == 2) {
1455     if (tracklet->validProj(0)) {
1456       addLayerProj(tracklet, 1);
1457     }
1458   }
1459 
1460   for (unsigned int disk = 2; disk < 6; disk++) {
1461     if (layer == 2 && disk == 5)
1462       continue;
1463     if (tracklet->validProj(N_LAYER + disk - 1)) {
1464       addDiskProj(tracklet, disk);
1465     }
1466   }
1467 
1468   return true;
1469 }