Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:21:50

0001 #include "L1Trigger/TrackFindingTMTT/interface/Settings.h"
0002 #include "L1Trigger/TrackFindingTMTT/interface/DigitalTrack.h"
0003 #include "L1Trigger/TrackFindingTMTT/interface/L1fittedTrack.h"
0004 
0005 #include "DataFormats/Math/interface/deltaPhi.h"
0006 #include "FWCore/Utilities/interface/Exception.h"
0007 
0008 using namespace std;
0009 
0010 namespace tmtt {
0011 
0012   //=== Note configuration parameters.
0013 
0014   DigitalTrack::DigitalTrack(const Settings* settings, const string& fitterName, const L1fittedTrack* fitTrk)
0015       :
0016 
0017         // Digitization configuration parameters
0018         settings_(settings),
0019 
0020         // Number of phi sectors and phi nonants.
0021         numPhiSectors_(settings->numPhiSectors()),
0022         numPhiNonants_(settings->numPhiNonants()),
0023         // Phi sector and phi nonant width (radians)
0024         phiSectorWidth_(2. * M_PI / float(numPhiSectors_)),
0025         phiNonantWidth_(2. * M_PI / float(numPhiNonants_)),
0026         // Radius from beamline with respect to which stub r coord. is measured.
0027         chosenRofPhi_(settings->chosenRofPhi()),
0028 
0029         // Number of q/Pt bins in Hough  transform array.
0030         nbinsPt_((int)settings->houghNbinsPt()),
0031         invPtToDPhi_(settings->invPtToDphi()),
0032 
0033         // Info about fitted track
0034         fitterName_(fitterName),
0035         nHelixParams_(fitTrk->nHelixParam()),
0036 
0037         nlayers_(fitTrk->numLayers()),
0038         iPhiSec_(fitTrk->iPhiSec()),
0039         iEtaReg_(fitTrk->iEtaReg()),
0040         mBin_(int(fitTrk->cellLocationHT().first) - floor(settings_->houghNbinsPt() / 2)),
0041         cBin_(int(fitTrk->cellLocationHT().second) - floor(settings_->houghNbinsPhi() / 2)),
0042         mBinhelix_(int(fitTrk->cellLocationFit().first) - floor(settings_->houghNbinsPt() / 2)),
0043         cBinhelix_(int(fitTrk->cellLocationFit().second) - floor(settings_->houghNbinsPhi() / 2)),
0044         hitPattern_(fitTrk->hitPattern()),
0045         consistent_(fitTrk->consistentHTcell()),
0046         consistentSect_(fitTrk->consistentSector()),
0047         accepted_(fitTrk->accepted()),
0048 
0049         qOverPt_orig_(fitTrk->qOverPt()),
0050         oneOver2r_orig_(fitTrk->qOverPt() * invPtToDPhi_),
0051         d0_orig_(fitTrk->d0()),
0052         phi0_orig_(fitTrk->phi0()),
0053         tanLambda_orig_(fitTrk->tanLambda()),
0054         z0_orig_(fitTrk->z0()),
0055         chisquaredRphi_orig_(fitTrk->chi2rphi()),
0056         chisquaredRz_orig_(fitTrk->chi2rz()),
0057 
0058         // Same again with beam-spot constraint.
0059         qOverPt_bcon_orig_(fitTrk->qOverPt_bcon()),
0060         oneOver2r_bcon_orig_(fitTrk->qOverPt_bcon() * invPtToDPhi_),
0061         phi0_bcon_orig_(fitTrk->phi0_bcon()),
0062         chisquaredRphi_bcon_orig_(fitTrk->chi2rphi_bcon()) {
0063     // Get digitisation parameters for this particular track fitter.
0064     this->loadDigiCfg(fitterName);
0065 
0066     // Complete storage of info about fitted track & truth particle.
0067     double phiCentreSec0 = -M_PI / float(numPhiNonants_) + M_PI / float(numPhiSectors_);
0068     phiSectorCentre_ = phiSectorWidth_ * float(iPhiSec_) + phiCentreSec0;
0069     phi0rel_orig_ = reco::deltaPhi(fitTrk->phi0(), phiSectorCentre_);
0070     phi0rel_bcon_orig_ = reco::deltaPhi(fitTrk->phi0_bcon(), phiSectorCentre_);
0071 
0072     // FIX: Remove this BODGE once BCHI increased to 11 in KFstate.h
0073     if (chisquaredRphi_orig_ >= chisquaredRange_)
0074       chisquaredRphi_orig_ = chisquaredRange_ - 0.1;
0075     if (chisquaredRphi_bcon_orig_ >= chisquaredRange_)
0076       chisquaredRphi_bcon_orig_ = chisquaredRange_ - 0.1;
0077 
0078     // Associated truth, if any.
0079     const TP* tp = fitTrk->matchedTP();
0080     bool tpOK = (tp != nullptr);
0081     tp_tanLambda_ = tpOK ? tp->tanLambda() : 0;
0082     tp_qoverpt_ = tpOK ? tp->qOverPt() : 0;
0083     tp_pt_ = tpOK ? tp->pt() : 0;
0084     tp_d0_ = tpOK ? tp->d0() : 0;
0085     tp_eta_ = tpOK ? tp->eta() : 0;
0086     tp_phi0_ = tpOK ? tp->phi0() : 0;
0087     tp_z0_ = tpOK ? tp->z0() : 0;
0088     tp_index_ = tpOK ? tp->index() : -1;
0089     tp_useForAlgEff_ = tpOK ? tp->useForAlgEff() : false;
0090     tp_useForEff_ = tpOK ? tp->useForEff() : false;
0091     tp_pdgId_ = tpOK ? tp->pdgId() : 0;
0092 
0093     // Digitize track.
0094     this->makeDigitalTrack();
0095   }
0096 
0097   //=== Load digitisation configuration parameters for the specific track fitter being used here.
0098 
0099   void DigitalTrack::loadDigiCfg(const string& fitterName) {
0100     if (fitterName == "SimpleLR4") {
0101       // SimpleLR4 track fitter
0102       skipTrackDigi_ = settings_->slr_skipTrackDigi();
0103       oneOver2rBits_ = settings_->slr_oneOver2rBits();
0104       oneOver2rRange_ = settings_->slr_oneOver2rRange();
0105       d0Bits_ = settings_->slr_d0Bits();
0106       d0Range_ = settings_->slr_d0Range();
0107       phi0Bits_ = settings_->slr_phi0Bits();
0108       phi0Range_ = settings_->slr_phi0Range();
0109       z0Bits_ = settings_->slr_z0Bits();
0110       z0Range_ = settings_->slr_z0Range();
0111       tanLambdaBits_ = settings_->slr_tanlambdaBits();
0112       tanLambdaRange_ = settings_->slr_tanlambdaRange();
0113       chisquaredBits_ = settings_->slr_chisquaredBits();
0114       chisquaredRange_ = settings_->slr_chisquaredRange();
0115     } else {
0116       // KF track fitter
0117       // Also used for all other fitters, though unlikely to be correct them them ...
0118       if (fitterName == "KF4ParamsComb" || fitterName == "KF5ParamsComb" || fitterName == "KF4ParamsCombHLS") {
0119         skipTrackDigi_ = settings_->kf_skipTrackDigi();
0120       } else {
0121         skipTrackDigi_ = settings_->other_skipTrackDigi();  // Allows to skip digitisation for other fitters
0122       }
0123       oneOver2rBits_ = settings_->kf_oneOver2rBits();
0124       oneOver2rRange_ = settings_->kf_oneOver2rRange();
0125       d0Bits_ = settings_->kf_d0Bits();
0126       d0Range_ = settings_->kf_d0Range();
0127       phi0Bits_ = settings_->kf_phi0Bits();
0128       phi0Range_ = settings_->kf_phi0Range();
0129       z0Bits_ = settings_->kf_z0Bits();
0130       z0Range_ = settings_->kf_z0Range();
0131       tanLambdaBits_ = settings_->kf_tanlambdaBits();
0132       tanLambdaRange_ = settings_->kf_tanlambdaRange();
0133       chisquaredBits_ = settings_->kf_chisquaredBits();
0134       chisquaredRange_ = settings_->kf_chisquaredRange();
0135     }
0136 
0137     // Calculate multipliers to digitize the floating point numbers.
0138     oneOver2rMult_ = pow(2., oneOver2rBits_) / oneOver2rRange_;
0139     d0Mult_ = pow(2., d0Bits_) / d0Range_;
0140     phi0Mult_ = pow(2., phi0Bits_) / phi0Range_;
0141     z0Mult_ = pow(2., z0Bits_) / z0Range_;
0142     tanLambdaMult_ = pow(2., tanLambdaBits_) / tanLambdaRange_;
0143     chisquaredMult_ = pow(2., chisquaredBits_) / chisquaredRange_;
0144   }
0145 
0146   //=== Digitize track
0147 
0148   void DigitalTrack::makeDigitalTrack() {
0149     if (skipTrackDigi_) {
0150       // Optionally skip track digitisaton if done internally inside track fitting code, so
0151       // retain original helix params.
0152       iDigi_oneOver2r_ = 0;
0153       iDigi_d0_ = 0;
0154       iDigi_phi0rel_ = 0;
0155       iDigi_tanLambda_ = 0;
0156       iDigi_z0_ = 0;
0157       iDigi_chisquaredRphi_ = 0;
0158       iDigi_chisquaredRz_ = 0;
0159 
0160       iDigi_oneOver2r_bcon_ = 0;
0161       iDigi_phi0rel_bcon_ = 0;
0162       iDigi_chisquaredRphi_bcon_ = 0;
0163 
0164       oneOver2r_ = oneOver2r_orig_;
0165       qOverPt_ = qOverPt_orig_;
0166       d0_ = d0_orig_;
0167       phi0rel_ = phi0rel_orig_;
0168       phi0_ = phi0_orig_;
0169       tanLambda_ = tanLambda_orig_;
0170       z0_ = z0_orig_;
0171       chisquaredRphi_ = chisquaredRphi_orig_;
0172       chisquaredRz_ = chisquaredRz_orig_;
0173 
0174       // Same again with beam-spot constraint.
0175       oneOver2r_bcon_ = oneOver2r_bcon_orig_;
0176       qOverPt_bcon_ = qOverPt_bcon_orig_;
0177       phi0rel_bcon_ = phi0rel_bcon_orig_;
0178       phi0_bcon_ = phi0_bcon_orig_;
0179       chisquaredRphi_bcon_ = chisquaredRphi_bcon_orig_;
0180 
0181     } else {
0182       //--- Digitize variables
0183 
0184       iDigi_oneOver2r_ = floor(oneOver2r_orig_ * oneOver2rMult_);
0185       iDigi_d0_ = floor(d0_orig_ * d0Mult_);
0186       iDigi_phi0rel_ = floor(phi0rel_orig_ * phi0Mult_);
0187       iDigi_tanLambda_ = floor(tanLambda_orig_ * tanLambdaMult_);
0188       iDigi_z0_ = floor(z0_orig_ * z0Mult_);
0189       iDigi_chisquaredRphi_ = floor(chisquaredRphi_orig_ * chisquaredMult_);
0190       iDigi_chisquaredRz_ = floor(chisquaredRz_orig_ * chisquaredMult_);
0191 
0192       // If fitted declared track invalid, it will have set its chi2 to very large number.
0193       // So truncate it at maximum allowed by digitisation range.
0194       if (!accepted_) {
0195         iDigi_chisquaredRphi_ = pow(2., chisquaredBits_) - 1;
0196         iDigi_chisquaredRz_ = pow(2., chisquaredBits_) - 1;
0197       }
0198 
0199       // Same again with beam-spot constraint.
0200       iDigi_oneOver2r_bcon_ = floor(oneOver2r_bcon_orig_ * oneOver2rMult_);
0201       iDigi_phi0rel_bcon_ = floor(phi0rel_bcon_orig_ * phi0Mult_);
0202       iDigi_chisquaredRphi_bcon_ = floor(chisquaredRphi_bcon_orig_ * chisquaredMult_);
0203       if (!accepted_)
0204         iDigi_chisquaredRphi_bcon_ = pow(2., chisquaredBits_) - 1;
0205 
0206       //--- Determine floating point track params from digitized numbers (so with degraded resolution).
0207 
0208       oneOver2r_ = (iDigi_oneOver2r_ + 0.5) / oneOver2rMult_;
0209       qOverPt_ = oneOver2r_ / invPtToDPhi_;
0210       if (nHelixParams_ == 5) {
0211         d0_ = (iDigi_d0_ + 0.5) / d0Mult_;
0212       } else {
0213         d0_ = 0.;
0214       }
0215       phi0rel_ = (iDigi_phi0rel_ + 0.5) / phi0Mult_;
0216       phi0_ = reco::deltaPhi(phi0rel_, -phiSectorCentre_);
0217       tanLambda_ = (iDigi_tanLambda_ + 0.5) / tanLambdaMult_;
0218       z0_ = (iDigi_z0_ + 0.5) / z0Mult_;
0219       chisquaredRphi_ = (iDigi_chisquaredRphi_ + 0.5) / chisquaredMult_;
0220       chisquaredRz_ = (iDigi_chisquaredRz_ + 0.5) / chisquaredMult_;
0221 
0222       // Same again with beam-spot constraint.
0223       if (nHelixParams_ == 5) {
0224         oneOver2r_bcon_ = (iDigi_oneOver2r_bcon_ + 0.5) / oneOver2rMult_;
0225         qOverPt_bcon_ = oneOver2r_bcon_ / invPtToDPhi_;
0226         phi0rel_bcon_ = (iDigi_phi0rel_bcon_ + 0.5) / phi0Mult_;
0227         phi0_bcon_ = reco::deltaPhi(phi0rel_bcon_, -phiSectorCentre_);
0228         chisquaredRphi_bcon_ = (iDigi_chisquaredRphi_bcon_ + 0.5) / chisquaredMult_;
0229       } else {
0230         oneOver2r_bcon_ = oneOver2r_;
0231         qOverPt_bcon_ = qOverPt_;
0232         phi0rel_bcon_ = phi0rel_;
0233         phi0_bcon_ = phi0_;
0234         chisquaredRphi_bcon_ = chisquaredRphi_;
0235       }
0236 
0237       // Check that track coords. are within assumed digitization range.
0238       this->checkInRange();
0239 
0240       // Check that digitization followed by undigitization doesn't change results too much.
0241       this->checkAccuracy();
0242     }
0243   }
0244 
0245   //=== Check that stub coords. are within assumed digitization range.
0246 
0247   void DigitalTrack::checkInRange() const {
0248     if (accepted_) {  // Don't bother apply to tracks rejected by the fitter.
0249       if (std::abs(oneOver2r_orig_) >= 0.5 * oneOver2rRange_)
0250         throw cms::Exception("BadConfig")
0251             << "DigitalTrack: Track oneOver2r is out of assumed digitization range."
0252             << " |oneOver2r| = " << std::abs(oneOver2r_orig_) << " > " << 0.5 * oneOver2rRange_
0253             << "; Fitter=" << fitterName_ << "; track accepted = " << accepted_;
0254       if (consistentSect_) {  // don't bother if track will fail sector consistency cut.
0255         if (std::abs(phi0rel_orig_) >= 0.5 * phi0Range_)
0256           throw cms::Exception("BadConfig") << "DigitalTrack: Track phi0rel is out of assumed digitization range."
0257                                             << " |phi0rel| = " << std::abs(phi0rel_orig_) << " > " << 0.5 * phi0Range_
0258                                             << "; Fitter=" << fitterName_ << "; track accepted = " << accepted_;
0259       }
0260       if (std::abs(z0_orig_) >= 0.5 * z0Range_)
0261         throw cms::Exception("BadConfig") << "DigitalTrack:  Track z0 is out of assumed digitization range."
0262                                           << " |z0| = " << std::abs(z0_orig_) << " > " << 0.5 * z0Range_
0263                                           << "; Fitter=" << fitterName_ << "; track accepted = " << accepted_;
0264       if (std::abs(d0_orig_) >= 0.5 * d0Range_)
0265         throw cms::Exception("BadConfig") << "DigitalTrack:  Track d0 is out of assumed digitization range."
0266                                           << " |d0| = " << std::abs(d0_orig_) << " > " << 0.5 * d0Range_
0267                                           << "; Fitter=" << fitterName_ << "; track accepted = " << accepted_;
0268       if (std::abs(tanLambda_orig_) >= 0.5 * tanLambdaRange_)
0269         throw cms::Exception("BadConfig")
0270             << "DigitalTrack: Track tanLambda is out of assumed digitization range."
0271             << " |tanLambda| = " << std::abs(tanLambda_orig_) << " > " << 0.5 * tanLambdaRange_
0272             << "; Fitter=" << fitterName_ << "; track accepted = " << accepted_;
0273       if (accepted_) {  // Tracks declared invalid by fitter can have very large original chi2.
0274         if (chisquaredRphi_orig_ >= chisquaredRange_ or chisquaredRphi_orig_ < 0.)
0275           throw cms::Exception("BadConfig")
0276               << "DigitalTrack: Track chisquaredRphi is out of assumed digitization range."
0277               << " chisquaredRphi = " << chisquaredRphi_orig_ << " > " << chisquaredRange_ << " or < 0"
0278               << "; Fitter=" << fitterName_ << "; track accepted = " << accepted_;
0279         if (chisquaredRz_orig_ >= chisquaredRange_ or chisquaredRz_orig_ < 0.)
0280           throw cms::Exception("BadConfig")
0281               << "DigitalTrack: Track chisquaredRz is out of assumed digitization range."
0282               << " chisquaredRz = " << chisquaredRz_orig_ << " > " << chisquaredRange_ << " or < 0"
0283               << "; Fitter=" << fitterName_ << "; track accepted = " << accepted_;
0284       }
0285     }
0286   }
0287 
0288   //=== Check that digitisation followed by undigitisation doesn't change significantly the stub coordinates.
0289 
0290   void DigitalTrack::checkAccuracy() const {
0291     if (accepted_) {  // Don't bother apply to tracks rejected by the fitter.
0292       float TA = qOverPt_ - qOverPt_orig_;
0293       float TB = reco::deltaPhi(phi0_, phi0_orig_);
0294       float TC = z0_ - z0_orig_;
0295       float TD = tanLambda_ - tanLambda_orig_;
0296       float TE = d0_ - d0_orig_;
0297       float TF = chisquaredRphi_ - chisquaredRphi_orig_;
0298       float TG = chisquaredRz_ - chisquaredRz_orig_;
0299 
0300       // Compare to small numbers, representing acceptable precision loss.
0301       constexpr float smallTA = 0.01, smallTB = 0.001, smallTC = 0.05, smallTD = 0.002, smallTE = 0.05, smallTF = 0.5,
0302                       smallTG = 0.5;
0303       if (std::abs(TA) > smallTA || std::abs(TB) > smallTB || std::abs(TC) > smallTC || std::abs(TD) > smallTD ||
0304           std::abs(TE) > smallTE || std::abs(TF) > smallTF || std::abs(TG) > smallTG) {
0305         throw cms::Exception("LogicError")
0306             << "WARNING: DigitalTrack lost precision: " << fitterName_ << " accepted=" << accepted_ << " " << TA << " "
0307             << TB << " " << TC << " " << TD << " " << TE << " " << TF << " " << TG;
0308       }
0309     }
0310   }
0311 
0312 }  // namespace tmtt