Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-12-02 23:44:55

0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/stream/EDProducer.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/EventSetup.h"
0005 #include "FWCore/Framework/interface/ESHandle.h"
0006 #include "FWCore/Framework/interface/ConsumesCollector.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 
0009 #include "RecoMTD/DetLayers/interface/MTDDetLayerGeometry.h"
0010 #include "RecoMTD/Records/interface/MTDRecoGeometryRecord.h"
0011 
0012 #include "MagneticField/Engine/interface/MagneticField.h"
0013 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0014 
0015 #include "TrackingTools/KalmanUpdators/interface/Chi2MeasurementEstimator.h"
0016 
0017 #include "DataFormats/TrackerRecHit2D/interface/MTDTrackingRecHit.h"
0018 
0019 #include "RecoMTD/DetLayers/interface/MTDTrayBarrelLayer.h"
0020 #include "TrackingTools/DetLayers/interface/ForwardDetLayer.h"
0021 
0022 #include "DataFormats/ForwardDetId/interface/BTLDetId.h"
0023 #include "DataFormats/ForwardDetId/interface/ETLDetId.h"
0024 #include "DataFormats/ForwardDetId/interface/MTDChannelIdentifier.h"
0025 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
0026 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0027 
0028 #include "TrackingTools/PatternTools/interface/Trajectory.h"
0029 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
0030 
0031 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0032 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0033 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0034 
0035 #include "TrackingTools/TransientTrackingRecHit/interface/TransientTrackingRecHit.h"
0036 
0037 #include "RecoMTD/TransientTrackingRecHit/interface/MTDTransientTrackingRecHitBuilder.h"
0038 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
0039 
0040 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
0041 #include "TrackingTools/GeomPropagators/interface/Propagator.h"
0042 
0043 #include "TrackingTools/PatternTools/interface/TSCBLBuilderWithPropagator.h"
0044 
0045 #include "RecoTracker/TransientTrackingRecHit/interface/Traj2TrackHits.h"
0046 #include "TrackingTools/TrackRefitter/interface/TrackTransformer.h"
0047 
0048 #include <sstream>
0049 
0050 #include "Geometry/CommonTopologies/interface/Topology.h"
0051 
0052 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0053 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0054 
0055 #include "DataFormats/Math/interface/GeantUnits.h"
0056 #include "DataFormats/Math/interface/LorentzVector.h"
0057 #include "CLHEP/Units/GlobalPhysicalConstants.h"
0058 #include "DataFormats/Math/interface/Rounding.h"
0059 
0060 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0061 #include "DataFormats/VertexReco/interface/Vertex.h"
0062 
0063 using namespace std;
0064 using namespace edm;
0065 using namespace reco;
0066 
0067 namespace {
0068   constexpr float c_cm_ns = geant_units::operators::convertMmToCm(CLHEP::c_light);  // [mm/ns] -> [cm/ns]
0069   constexpr float c_inv = 1.0f / c_cm_ns;
0070 
0071   class MTDHitMatchingInfo {
0072   public:
0073     MTDHitMatchingInfo() {
0074       hit = nullptr;
0075       estChi2 = std::numeric_limits<float>::max();
0076       timeChi2 = std::numeric_limits<float>::max();
0077     }
0078 
0079     //Operator used to sort the hits while performing the matching step at the MTD
0080     inline bool operator<(const MTDHitMatchingInfo& m2) const {
0081       //only for good matching in time use estChi2, otherwise use mostly time compatibility
0082       constexpr float chi2_cut = 10.f;
0083       constexpr float low_weight = 3.f;
0084       constexpr float high_weight = 8.f;
0085       if (timeChi2 < chi2_cut && m2.timeChi2 < chi2_cut)
0086         return chi2(low_weight) < m2.chi2(low_weight);
0087       else
0088         return chi2(high_weight) < m2.chi2(high_weight);
0089     }
0090 
0091     inline float chi2(float timeWeight = 1.f) const { return estChi2 + timeWeight * timeChi2; }
0092 
0093     const MTDTrackingRecHit* hit;
0094     float estChi2;
0095     float timeChi2;
0096   };
0097 
0098   class TrackSegments {
0099   public:
0100     TrackSegments() {
0101       sigmaTofs_.reserve(30);  // observed upper limit on nSegments
0102     };
0103 
0104     inline uint32_t addSegment(float tPath, float tMom2, float sigmaMom) {
0105       segmentPathOvc_.emplace_back(tPath * c_inv);
0106       segmentMom2_.emplace_back(tMom2);
0107       segmentSigmaMom_.emplace_back(sigmaMom);
0108       nSegment_++;
0109 
0110       LogTrace("TrackExtenderWithMTD") << "addSegment # " << nSegment_ << " s = " << tPath
0111                                        << " p = " << std::sqrt(tMom2) << " sigma_p = " << sigmaMom
0112                                        << " sigma_p/p = " << sigmaMom / std::sqrt(tMom2) * 100 << " %";
0113 
0114       return nSegment_;
0115     }
0116 
0117     inline float computeTof(float mass_inv2) const {
0118       float tof(0.f);
0119       for (uint32_t iSeg = 0; iSeg < nSegment_; iSeg++) {
0120         float gammasq = 1.f + segmentMom2_[iSeg] * mass_inv2;
0121         float beta = std::sqrt(1.f - 1.f / gammasq);
0122         tof += segmentPathOvc_[iSeg] / beta;
0123 
0124         LogTrace("TrackExtenderWithMTD") << " TOF Segment # " << iSeg + 1 << " p = " << std::sqrt(segmentMom2_[iSeg])
0125                                          << " tof = " << tof;
0126 
0127 #ifdef EDM_ML_DEBUG
0128         float sigma_tof = segmentPathOvc_[iSeg] * segmentSigmaMom_[iSeg] /
0129                           (segmentMom2_[iSeg] * sqrt(segmentMom2_[iSeg] + 1 / mass_inv2) * mass_inv2);
0130 
0131         LogTrace("TrackExtenderWithMTD") << "TOF Segment # " << iSeg + 1 << std::fixed << std::setw(6)
0132                                          << " tof segment = " << segmentPathOvc_[iSeg] / beta << std::scientific
0133                                          << "+/- " << sigma_tof << std::fixed
0134                                          << "(rel. err. = " << sigma_tof / (segmentPathOvc_[iSeg] / beta) * 100
0135                                          << " %)";
0136 #endif
0137       }
0138 
0139       return tof;
0140     }
0141 
0142     inline float computeSigmaTof(float mass_inv2) {
0143       float sigmatof = 0.;
0144 
0145       // remove previously calculated sigmaTofs
0146       sigmaTofs_.clear();
0147 
0148       // compute sigma(tof) on each segment first by propagating sigma(p)
0149       // also add diagonal terms to sigmatof
0150       float sigma = 0.;
0151       for (uint32_t iSeg = 0; iSeg < nSegment_; iSeg++) {
0152         sigma = segmentPathOvc_[iSeg] * segmentSigmaMom_[iSeg] /
0153                 (segmentMom2_[iSeg] * sqrt(segmentMom2_[iSeg] + 1 / mass_inv2) * mass_inv2);
0154         sigmaTofs_.push_back(sigma);
0155 
0156         sigmatof += sigma * sigma;
0157       }
0158 
0159       // compute sigma on sum of tofs assuming full correlation between segments
0160       for (uint32_t iSeg = 0; iSeg < nSegment_; iSeg++) {
0161         for (uint32_t jSeg = iSeg + 1; jSeg < nSegment_; jSeg++) {
0162           sigmatof += 2 * sigmaTofs_[iSeg] * sigmaTofs_[jSeg];
0163         }
0164       }
0165 
0166       return sqrt(sigmatof);
0167     }
0168 
0169     inline uint32_t size() const { return nSegment_; }
0170 
0171     inline uint32_t removeFirstSegment() {
0172       if (nSegment_ > 0) {
0173         segmentPathOvc_.erase(segmentPathOvc_.begin());
0174         segmentMom2_.erase(segmentMom2_.begin());
0175         nSegment_--;
0176       }
0177       return nSegment_;
0178     }
0179 
0180     inline std::pair<float, float> getSegmentPathAndMom2(uint32_t iSegment) const {
0181       if (iSegment >= nSegment_) {
0182         throw cms::Exception("TrackExtenderWithMTD") << "Requesting non existing track segment #" << iSegment;
0183       }
0184       return std::make_pair(segmentPathOvc_[iSegment], segmentMom2_[iSegment]);
0185     }
0186 
0187     uint32_t nSegment_ = 0;
0188     std::vector<float> segmentPathOvc_;
0189     std::vector<float> segmentMom2_;
0190     std::vector<float> segmentSigmaMom_;
0191 
0192     std::vector<float> sigmaTofs_;
0193   };
0194 
0195   struct TrackTofPidInfo {
0196     float tmtd;
0197     float tmtderror;
0198     float pathlength;
0199 
0200     float betaerror;
0201 
0202     float dt;
0203     float dterror;
0204     float dterror2;
0205     float dtchi2;
0206 
0207     float dt_best;
0208     float dterror_best;
0209     float dtchi2_best;
0210 
0211     float gammasq_pi;
0212     float beta_pi;
0213     float dt_pi;
0214     float sigma_dt_pi;
0215 
0216     float gammasq_k;
0217     float beta_k;
0218     float dt_k;
0219     float sigma_dt_k;
0220 
0221     float gammasq_p;
0222     float beta_p;
0223     float dt_p;
0224     float sigma_dt_p;
0225 
0226     float prob_pi;
0227     float prob_k;
0228     float prob_p;
0229   };
0230 
0231   enum class TofCalc { kCost = 1, kSegm = 2, kMixd = 3 };
0232   enum class SigmaTofCalc { kCost = 1, kSegm = 2, kMixd = 3 };
0233 
0234   const TrackTofPidInfo computeTrackTofPidInfo(float magp2,
0235                                                float length,
0236                                                TrackSegments trs,
0237                                                float t_mtd,
0238                                                float t_mtderr,
0239                                                float t_vtx,
0240                                                float t_vtx_err,
0241                                                bool addPIDError = true,
0242                                                TofCalc choice = TofCalc::kCost,
0243                                                SigmaTofCalc sigma_choice = SigmaTofCalc::kCost) {
0244     constexpr float m_pi = 0.13957018f;
0245     constexpr float m_pi_inv2 = 1.0f / m_pi / m_pi;
0246     constexpr float m_k = 0.493677f;
0247     constexpr float m_k_inv2 = 1.0f / m_k / m_k;
0248     constexpr float m_p = 0.9382720813f;
0249     constexpr float m_p_inv2 = 1.0f / m_p / m_p;
0250 
0251     TrackTofPidInfo tofpid;
0252 
0253     tofpid.tmtd = t_mtd;
0254     tofpid.tmtderror = t_mtderr;
0255     tofpid.pathlength = length;
0256 
0257     auto deltat = [&](const float mass_inv2, const float betatmp) {
0258       float res(1.f);
0259       switch (choice) {
0260         case TofCalc::kCost:
0261           res = tofpid.pathlength / betatmp * c_inv;
0262           break;
0263         case TofCalc::kSegm:
0264           res = trs.computeTof(mass_inv2);
0265           break;
0266         case TofCalc::kMixd:
0267           res = trs.computeTof(mass_inv2) + tofpid.pathlength / betatmp * c_inv;
0268           break;
0269       }
0270       return res;
0271     };
0272 
0273     auto sigmadeltat = [&](const float mass_inv2) {
0274       float res(1.f);
0275       switch (sigma_choice) {
0276         case SigmaTofCalc::kCost:
0277           // sigma(t) = sigma(p) * |dt/dp| = sigma(p) * DeltaL/c * m^2 / (p^2 * E)
0278           res = tofpid.pathlength * c_inv * trs.segmentSigmaMom_[trs.nSegment_ - 1] /
0279                 (magp2 * sqrt(magp2 + 1 / mass_inv2) * mass_inv2);
0280           break;
0281         case SigmaTofCalc::kSegm:
0282           res = trs.computeSigmaTof(mass_inv2);
0283           break;
0284         case SigmaTofCalc::kMixd:
0285           float res1 = tofpid.pathlength * c_inv * trs.segmentSigmaMom_[trs.nSegment_ - 1] /
0286                        (magp2 * sqrt(magp2 + 1 / mass_inv2) * mass_inv2);
0287           float res2 = trs.computeSigmaTof(mass_inv2);
0288           res = sqrt(res1 * res1 + res2 * res2 + 2 * res1 * res2);
0289       }
0290 
0291       return res;
0292     };
0293 
0294     tofpid.gammasq_pi = 1.f + magp2 * m_pi_inv2;
0295     tofpid.beta_pi = std::sqrt(1.f - 1.f / tofpid.gammasq_pi);
0296     tofpid.dt_pi = deltat(m_pi_inv2, tofpid.beta_pi);
0297     tofpid.sigma_dt_pi = sigmadeltat(m_pi_inv2);
0298 
0299     tofpid.gammasq_k = 1.f + magp2 * m_k_inv2;
0300     tofpid.beta_k = std::sqrt(1.f - 1.f / tofpid.gammasq_k);
0301     tofpid.dt_k = deltat(m_k_inv2, tofpid.beta_k);
0302     tofpid.sigma_dt_k = sigmadeltat(m_k_inv2);
0303 
0304     tofpid.gammasq_p = 1.f + magp2 * m_p_inv2;
0305     tofpid.beta_p = std::sqrt(1.f - 1.f / tofpid.gammasq_p);
0306     tofpid.dt_p = deltat(m_p_inv2, tofpid.beta_p);
0307     tofpid.sigma_dt_p = sigmadeltat(m_p_inv2);
0308 
0309     tofpid.dt = tofpid.tmtd - tofpid.dt_pi - t_vtx;  //assume by default the pi hypothesis
0310     tofpid.dterror2 = tofpid.tmtderror * tofpid.tmtderror + t_vtx_err * t_vtx_err;
0311     tofpid.betaerror = 0.f;
0312     if (addPIDError) {
0313       tofpid.dterror2 = tofpid.dterror2 + (tofpid.dt_p - tofpid.dt_pi) * (tofpid.dt_p - tofpid.dt_pi);
0314       tofpid.betaerror = tofpid.beta_p - tofpid.beta_pi;
0315     } else {
0316       // only add sigma(TOF) if not considering mass hp. uncertainty
0317       tofpid.dterror2 = tofpid.dterror2 + tofpid.sigma_dt_pi * tofpid.sigma_dt_pi;
0318     }
0319     tofpid.dterror = sqrt(tofpid.dterror2);
0320 
0321     tofpid.dtchi2 = (tofpid.dt * tofpid.dt) / tofpid.dterror2;
0322 
0323     tofpid.dt_best = tofpid.dt;
0324     tofpid.dterror_best = tofpid.dterror;
0325     tofpid.dtchi2_best = tofpid.dtchi2;
0326 
0327     tofpid.prob_pi = -1.f;
0328     tofpid.prob_k = -1.f;
0329     tofpid.prob_p = -1.f;
0330 
0331     if (!addPIDError) {
0332       //*TODO* deal with heavier nucleons and/or BSM case here?
0333       const float dterror2_wo_sigmatof = tofpid.dterror2 - tofpid.sigma_dt_pi * tofpid.sigma_dt_pi;
0334       float chi2_pi = tofpid.dtchi2;
0335       float chi2_k = (tofpid.tmtd - tofpid.dt_k - t_vtx) * (tofpid.tmtd - tofpid.dt_k - t_vtx) /
0336                      (dterror2_wo_sigmatof + tofpid.sigma_dt_k * tofpid.sigma_dt_k);
0337       float chi2_p = (tofpid.tmtd - tofpid.dt_p - t_vtx) * (tofpid.tmtd - tofpid.dt_p - t_vtx) /
0338                      (dterror2_wo_sigmatof + tofpid.sigma_dt_p * tofpid.sigma_dt_p);
0339 
0340       float rawprob_pi = exp(-0.5f * chi2_pi);
0341       float rawprob_k = exp(-0.5f * chi2_k);
0342       float rawprob_p = exp(-0.5f * chi2_p);
0343       float normprob = 1.f / (rawprob_pi + rawprob_k + rawprob_p);
0344 
0345       tofpid.prob_pi = rawprob_pi * normprob;
0346       tofpid.prob_k = rawprob_k * normprob;
0347       tofpid.prob_p = rawprob_p * normprob;
0348 
0349       float prob_heavy = 1.f - tofpid.prob_pi;
0350       constexpr float heavy_threshold = 0.75f;
0351 
0352       if (prob_heavy > heavy_threshold) {
0353         if (chi2_k < chi2_p) {
0354           tofpid.dt_best = (tofpid.tmtd - tofpid.dt_k - t_vtx);
0355           tofpid.dtchi2_best = chi2_k;
0356         } else {
0357           tofpid.dt_best = (tofpid.tmtd - tofpid.dt_p - t_vtx);
0358           tofpid.dtchi2_best = chi2_p;
0359         }
0360       }
0361     }
0362     return tofpid;
0363   }
0364 
0365   bool getTrajectoryStateClosestToBeamLine(const Trajectory& traj,
0366                                            const reco::BeamSpot& bs,
0367                                            const Propagator* thePropagator,
0368                                            TrajectoryStateClosestToBeamLine& tscbl) {
0369     // get the state closest to the beamline
0370     TrajectoryStateOnSurface stateForProjectionToBeamLineOnSurface =
0371         traj.closestMeasurement(GlobalPoint(bs.x0(), bs.y0(), bs.z0())).updatedState();
0372 
0373     if (!stateForProjectionToBeamLineOnSurface.isValid()) {
0374       edm::LogError("CannotPropagateToBeamLine") << "the state on the closest measurement isnot valid. skipping track.";
0375       return false;
0376     }
0377 
0378     const FreeTrajectoryState& stateForProjectionToBeamLine = *stateForProjectionToBeamLineOnSurface.freeState();
0379 
0380     TSCBLBuilderWithPropagator tscblBuilder(*thePropagator);
0381     tscbl = tscblBuilder(stateForProjectionToBeamLine, bs);
0382 
0383     return tscbl.isValid();
0384   }
0385 
0386   bool trackPathLength(const Trajectory& traj,
0387                        const TrajectoryStateClosestToBeamLine& tscbl,
0388                        const Propagator* thePropagator,
0389                        float& pathlength,
0390                        TrackSegments& trs) {
0391     pathlength = 0.f;
0392 
0393     bool validpropagation = true;
0394     float oldp = traj.measurements().begin()->updatedState().globalMomentum().mag();
0395     float pathlength1 = 0.f;
0396     float pathlength2 = 0.f;
0397 
0398     //add pathlength layer by layer
0399     for (auto it = traj.measurements().begin(); it != traj.measurements().end() - 1; ++it) {
0400       const auto& propresult = thePropagator->propagateWithPath(it->updatedState(), (it + 1)->updatedState().surface());
0401       float layerpathlength = std::abs(propresult.second);
0402       if (layerpathlength == 0.f) {
0403         validpropagation = false;
0404       }
0405       pathlength1 += layerpathlength;
0406 
0407       // sigma(p) from curvilinear error (on q/p)
0408       float sigma_p = sqrt((it + 1)->updatedState().curvilinearError().matrix()(0, 0)) *
0409                       (it + 1)->updatedState().globalMomentum().mag2();
0410 
0411       trs.addSegment(layerpathlength, (it + 1)->updatedState().globalMomentum().mag2(), sigma_p);
0412 
0413       LogTrace("TrackExtenderWithMTD") << "TSOS " << std::fixed << std::setw(4) << trs.size() << " R_i " << std::fixed
0414                                        << std::setw(14) << it->updatedState().globalPosition().perp() << " z_i "
0415                                        << std::fixed << std::setw(14) << it->updatedState().globalPosition().z()
0416                                        << " R_e " << std::fixed << std::setw(14)
0417                                        << (it + 1)->updatedState().globalPosition().perp() << " z_e " << std::fixed
0418                                        << std::setw(14) << (it + 1)->updatedState().globalPosition().z() << " p "
0419                                        << std::fixed << std::setw(14) << (it + 1)->updatedState().globalMomentum().mag()
0420                                        << " dp " << std::fixed << std::setw(14)
0421                                        << (it + 1)->updatedState().globalMomentum().mag() - oldp;
0422       oldp = (it + 1)->updatedState().globalMomentum().mag();
0423     }
0424 
0425     //add distance from bs to first measurement
0426     auto const& tscblPCA = tscbl.trackStateAtPCA();
0427     auto const& aSurface = traj.direction() == alongMomentum ? traj.firstMeasurement().updatedState().surface()
0428                                                              : traj.lastMeasurement().updatedState().surface();
0429     pathlength2 = thePropagator->propagateWithPath(tscblPCA, aSurface).second;
0430     if (pathlength2 == 0.f) {
0431       validpropagation = false;
0432     }
0433     pathlength = pathlength1 + pathlength2;
0434 
0435     float sigma_p = sqrt(tscblPCA.curvilinearError().matrix()(0, 0)) * tscblPCA.momentum().mag2();
0436 
0437     trs.addSegment(pathlength2, tscblPCA.momentum().mag2(), sigma_p);
0438 
0439     LogTrace("TrackExtenderWithMTD") << "TSOS " << std::fixed << std::setw(4) << trs.size() << " R_e " << std::fixed
0440                                      << std::setw(14) << tscblPCA.position().perp() << " z_e " << std::fixed
0441                                      << std::setw(14) << tscblPCA.position().z() << " p " << std::fixed << std::setw(14)
0442                                      << tscblPCA.momentum().mag() << " dp " << std::fixed << std::setw(14)
0443                                      << tscblPCA.momentum().mag() - oldp << " sigma_p = " << std::fixed << std::setw(14)
0444                                      << sigma_p << " sigma_p/p = " << std::fixed << std::setw(14)
0445                                      << sigma_p / tscblPCA.momentum().mag() * 100 << " %";
0446 
0447     return validpropagation;
0448   }
0449 
0450   bool trackPathLength(const Trajectory& traj,
0451                        const reco::BeamSpot& bs,
0452                        const Propagator* thePropagator,
0453                        float& pathlength,
0454                        TrackSegments& trs) {
0455     pathlength = 0.f;
0456 
0457     TrajectoryStateClosestToBeamLine tscbl;
0458     bool tscbl_status = getTrajectoryStateClosestToBeamLine(traj, bs, thePropagator, tscbl);
0459 
0460     if (!tscbl_status)
0461       return false;
0462 
0463     return trackPathLength(traj, tscbl, thePropagator, pathlength, trs);
0464   }
0465 
0466 }  // namespace
0467 
0468 template <class TrackCollection>
0469 class TrackExtenderWithMTDT : public edm::stream::EDProducer<> {
0470 public:
0471   typedef typename TrackCollection::value_type TrackType;
0472   typedef edm::View<TrackType> InputCollection;
0473 
0474   TrackExtenderWithMTDT(const ParameterSet& pset);
0475 
0476   template <class H, class T>
0477   void fillValueMap(edm::Event& iEvent, const H& handle, const std::vector<T>& vec, const edm::EDPutToken& token) const;
0478 
0479   void produce(edm::Event& ev, const edm::EventSetup& es) final;
0480 
0481   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0482 
0483   TransientTrackingRecHit::ConstRecHitContainer tryBTLLayers(const TrajectoryStateOnSurface&,
0484                                                              const Trajectory& traj,
0485                                                              const float,
0486                                                              const float,
0487                                                              const TrackSegments&,
0488                                                              const MTDTrackingDetSetVector&,
0489                                                              const MTDDetLayerGeometry*,
0490                                                              const MagneticField* field,
0491                                                              const Propagator* prop,
0492                                                              const reco::BeamSpot& bs,
0493                                                              const float vtxTime,
0494                                                              const bool matchVertex,
0495                                                              MTDHitMatchingInfo& bestHit) const;
0496 
0497   TransientTrackingRecHit::ConstRecHitContainer tryETLLayers(const TrajectoryStateOnSurface&,
0498                                                              const Trajectory& traj,
0499                                                              const float,
0500                                                              const float,
0501                                                              const TrackSegments&,
0502                                                              const MTDTrackingDetSetVector&,
0503                                                              const MTDDetLayerGeometry*,
0504                                                              const MagneticField* field,
0505                                                              const Propagator* prop,
0506                                                              const reco::BeamSpot& bs,
0507                                                              const float vtxTime,
0508                                                              const bool matchVertex,
0509                                                              MTDHitMatchingInfo& bestHit) const;
0510 
0511   void fillMatchingHits(const DetLayer*,
0512                         const TrajectoryStateOnSurface&,
0513                         const Trajectory&,
0514                         const float,
0515                         const float,
0516                         const TrackSegments&,
0517                         const MTDTrackingDetSetVector&,
0518                         const Propagator*,
0519                         const reco::BeamSpot&,
0520                         const float&,
0521                         const bool,
0522                         TransientTrackingRecHit::ConstRecHitContainer&,
0523                         MTDHitMatchingInfo&) const;
0524 
0525   RefitDirection::GeometricalDirection checkRecHitsOrdering(
0526       TransientTrackingRecHit::ConstRecHitContainer const& recHits) const {
0527     if (!recHits.empty()) {
0528       GlobalPoint first = gtg_->idToDet(recHits.front()->geographicalId())->position();
0529       GlobalPoint last = gtg_->idToDet(recHits.back()->geographicalId())->position();
0530 
0531       // maybe perp2?
0532       auto rFirst = first.mag2();
0533       auto rLast = last.mag2();
0534       if (rFirst < rLast)
0535         return RefitDirection::insideOut;
0536       if (rFirst > rLast)
0537         return RefitDirection::outsideIn;
0538     }
0539     LogDebug("TrackExtenderWithMTD") << "Impossible to determine the rechits order" << endl;
0540     return RefitDirection::undetermined;
0541   }
0542 
0543   reco::Track buildTrack(const reco::TrackRef&,
0544                          const Trajectory&,
0545                          const Trajectory&,
0546                          const reco::BeamSpot&,
0547                          const MagneticField* field,
0548                          const Propagator* prop,
0549                          bool hasMTD,
0550                          float& pathLength,
0551                          float& tmtdOut,
0552                          float& sigmatmtdOut,
0553                          GlobalPoint& tmtdPosOut,
0554                          float& tofpi,
0555                          float& tofk,
0556                          float& tofp,
0557                          float& sigmatofpi,
0558                          float& sigmatofk,
0559                          float& sigmatofp) const;
0560   reco::TrackExtra buildTrackExtra(const Trajectory& trajectory) const;
0561 
0562   string dumpLayer(const DetLayer* layer) const;
0563 
0564 private:
0565   edm::EDPutToken btlMatchChi2Token_;
0566   edm::EDPutToken etlMatchChi2Token_;
0567   edm::EDPutToken btlMatchTimeChi2Token_;
0568   edm::EDPutToken etlMatchTimeChi2Token_;
0569   edm::EDPutToken npixBarrelToken_;
0570   edm::EDPutToken npixEndcapToken_;
0571   edm::EDPutToken outermostHitPositionToken_;
0572   edm::EDPutToken pOrigTrkToken_;
0573   edm::EDPutToken betaOrigTrkToken_;
0574   edm::EDPutToken t0OrigTrkToken_;
0575   edm::EDPutToken sigmat0OrigTrkToken_;
0576   edm::EDPutToken pathLengthOrigTrkToken_;
0577   edm::EDPutToken tmtdOrigTrkToken_;
0578   edm::EDPutToken sigmatmtdOrigTrkToken_;
0579   edm::EDPutToken tmtdPosOrigTrkToken_;
0580   edm::EDPutToken tofpiOrigTrkToken_;
0581   edm::EDPutToken tofkOrigTrkToken_;
0582   edm::EDPutToken tofpOrigTrkToken_;
0583   edm::EDPutToken sigmatofpiOrigTrkToken_;
0584   edm::EDPutToken sigmatofkOrigTrkToken_;
0585   edm::EDPutToken sigmatofpOrigTrkToken_;
0586   edm::EDPutToken assocOrigTrkToken_;
0587 
0588   edm::EDGetTokenT<InputCollection> tracksToken_;
0589   edm::EDGetTokenT<TrajTrackAssociationCollection> trajTrackAToken_;
0590   edm::EDGetTokenT<MTDTrackingDetSetVector> hitsToken_;
0591   edm::EDGetTokenT<reco::BeamSpot> bsToken_;
0592   edm::EDGetTokenT<GlobalPoint> genVtxPositionToken_;
0593   edm::EDGetTokenT<float> genVtxTimeToken_;
0594   edm::EDGetTokenT<VertexCollection> vtxToken_;
0595 
0596   const bool updateTraj_, updateExtra_, updatePattern_;
0597   const std::string mtdRecHitBuilder_, propagator_, transientTrackBuilder_;
0598   std::unique_ptr<MeasurementEstimator> theEstimator;
0599   std::unique_ptr<TrackTransformer> theTransformer;
0600   edm::ESHandle<TransientTrackBuilder> builder_;
0601   edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> builderToken_;
0602   edm::ESHandle<TransientTrackingRecHitBuilder> hitbuilder_;
0603   edm::ESGetToken<TransientTrackingRecHitBuilder, TransientRecHitRecord> hitbuilderToken_;
0604   edm::ESHandle<GlobalTrackingGeometry> gtg_;
0605   edm::ESGetToken<GlobalTrackingGeometry, GlobalTrackingGeometryRecord> gtgToken_;
0606 
0607   edm::ESGetToken<MTDDetLayerGeometry, MTDRecoGeometryRecord> dlgeoToken_;
0608   edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> magfldToken_;
0609   edm::ESGetToken<Propagator, TrackingComponentsRecord> propToken_;
0610   edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> ttopoToken_;
0611 
0612   const float estMaxChi2_;
0613   const float estMaxNSigma_;
0614   const float btlChi2Cut_;
0615   const float btlTimeChi2Cut_;
0616   const float etlChi2Cut_;
0617   const float etlTimeChi2Cut_;
0618 
0619   const bool useVertex_;
0620   const bool useSimVertex_;
0621   const float dzCut_;
0622   const float bsTimeSpread_;
0623 
0624   static constexpr float trackMaxBtlEta_ = 1.5;
0625 };
0626 
0627 template <class TrackCollection>
0628 TrackExtenderWithMTDT<TrackCollection>::TrackExtenderWithMTDT(const ParameterSet& iConfig)
0629     : tracksToken_(consumes<InputCollection>(iConfig.getParameter<edm::InputTag>("tracksSrc"))),
0630       trajTrackAToken_(consumes<TrajTrackAssociationCollection>(iConfig.getParameter<edm::InputTag>("trjtrkAssSrc"))),
0631       hitsToken_(consumes<MTDTrackingDetSetVector>(iConfig.getParameter<edm::InputTag>("hitsSrc"))),
0632       bsToken_(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpotSrc"))),
0633       updateTraj_(iConfig.getParameter<bool>("updateTrackTrajectory")),
0634       updateExtra_(iConfig.getParameter<bool>("updateTrackExtra")),
0635       updatePattern_(iConfig.getParameter<bool>("updateTrackHitPattern")),
0636       mtdRecHitBuilder_(iConfig.getParameter<std::string>("MTDRecHitBuilder")),
0637       propagator_(iConfig.getParameter<std::string>("Propagator")),
0638       transientTrackBuilder_(iConfig.getParameter<std::string>("TransientTrackBuilder")),
0639       estMaxChi2_(iConfig.getParameter<double>("estimatorMaxChi2")),
0640       estMaxNSigma_(iConfig.getParameter<double>("estimatorMaxNSigma")),
0641       btlChi2Cut_(iConfig.getParameter<double>("btlChi2Cut")),
0642       btlTimeChi2Cut_(iConfig.getParameter<double>("btlTimeChi2Cut")),
0643       etlChi2Cut_(iConfig.getParameter<double>("etlChi2Cut")),
0644       etlTimeChi2Cut_(iConfig.getParameter<double>("etlTimeChi2Cut")),
0645       useVertex_(iConfig.getParameter<bool>("useVertex")),
0646       useSimVertex_(iConfig.getParameter<bool>("useSimVertex")),
0647       dzCut_(iConfig.getParameter<double>("dZCut")),
0648       bsTimeSpread_(iConfig.getParameter<double>("bsTimeSpread")) {
0649   if (useVertex_) {
0650     if (useSimVertex_) {
0651       genVtxPositionToken_ = consumes<GlobalPoint>(iConfig.getParameter<edm::InputTag>("genVtxPositionSrc"));
0652       genVtxTimeToken_ = consumes<float>(iConfig.getParameter<edm::InputTag>("genVtxTimeSrc"));
0653     } else
0654       vtxToken_ = consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vtxSrc"));
0655   }
0656 
0657   theEstimator = std::make_unique<Chi2MeasurementEstimator>(estMaxChi2_, estMaxNSigma_);
0658   theTransformer = std::make_unique<TrackTransformer>(iConfig.getParameterSet("TrackTransformer"), consumesCollector());
0659 
0660   btlMatchChi2Token_ = produces<edm::ValueMap<float>>("btlMatchChi2");
0661   etlMatchChi2Token_ = produces<edm::ValueMap<float>>("etlMatchChi2");
0662   btlMatchTimeChi2Token_ = produces<edm::ValueMap<float>>("btlMatchTimeChi2");
0663   etlMatchTimeChi2Token_ = produces<edm::ValueMap<float>>("etlMatchTimeChi2");
0664   npixBarrelToken_ = produces<edm::ValueMap<int>>("npixBarrel");
0665   npixEndcapToken_ = produces<edm::ValueMap<int>>("npixEndcap");
0666   outermostHitPositionToken_ = produces<edm::ValueMap<float>>("generalTrackOutermostHitPosition");
0667   pOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackp");
0668   betaOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackBeta");
0669   t0OrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackt0");
0670   sigmat0OrigTrkToken_ = produces<edm::ValueMap<float>>("generalTracksigmat0");
0671   pathLengthOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackPathLength");
0672   tmtdOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTracktmtd");
0673   sigmatmtdOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTracksigmatmtd");
0674   tmtdPosOrigTrkToken_ = produces<edm::ValueMap<GlobalPoint>>("generalTrackmtdpos");
0675   tofpiOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackTofPi");
0676   tofkOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackTofK");
0677   tofpOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackTofP");
0678   sigmatofpiOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackSigmaTofPi");
0679   sigmatofkOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackSigmaTofK");
0680   sigmatofpOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackSigmaTofP");
0681   assocOrigTrkToken_ = produces<edm::ValueMap<int>>("generalTrackassoc");
0682 
0683   builderToken_ = esConsumes<TransientTrackBuilder, TransientTrackRecord>(edm::ESInputTag("", transientTrackBuilder_));
0684   hitbuilderToken_ =
0685       esConsumes<TransientTrackingRecHitBuilder, TransientRecHitRecord>(edm::ESInputTag("", mtdRecHitBuilder_));
0686   gtgToken_ = esConsumes<GlobalTrackingGeometry, GlobalTrackingGeometryRecord>();
0687   dlgeoToken_ = esConsumes<MTDDetLayerGeometry, MTDRecoGeometryRecord>();
0688   magfldToken_ = esConsumes<MagneticField, IdealMagneticFieldRecord>();
0689   propToken_ = esConsumes<Propagator, TrackingComponentsRecord>(edm::ESInputTag("", propagator_));
0690   ttopoToken_ = esConsumes<TrackerTopology, TrackerTopologyRcd>();
0691 
0692   produces<edm::OwnVector<TrackingRecHit>>();
0693   produces<reco::TrackExtraCollection>();
0694   produces<TrackCollection>();
0695 }
0696 
0697 template <class TrackCollection>
0698 void TrackExtenderWithMTDT<TrackCollection>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0699   edm::ParameterSetDescription desc, transDesc;
0700   desc.add<edm::InputTag>("tracksSrc", edm::InputTag("generalTracks"));
0701   desc.add<edm::InputTag>("trjtrkAssSrc", edm::InputTag("generalTracks"));
0702   desc.add<edm::InputTag>("hitsSrc", edm::InputTag("mtdTrackingRecHits"));
0703   desc.add<edm::InputTag>("beamSpotSrc", edm::InputTag("offlineBeamSpot"));
0704   desc.add<edm::InputTag>("genVtxPositionSrc", edm::InputTag("genParticles:xyz0"));
0705   desc.add<edm::InputTag>("genVtxTimeSrc", edm::InputTag("genParticles:t0"));
0706   desc.add<edm::InputTag>("vtxSrc", edm::InputTag("offlinePrimaryVertices4D"));
0707   desc.add<bool>("updateTrackTrajectory", true);
0708   desc.add<bool>("updateTrackExtra", true);
0709   desc.add<bool>("updateTrackHitPattern", true);
0710   desc.add<std::string>("TransientTrackBuilder", "TransientTrackBuilder");
0711   desc.add<std::string>("MTDRecHitBuilder", "MTDRecHitBuilder");
0712   desc.add<std::string>("Propagator", "PropagatorWithMaterialForMTD");
0713   TrackTransformer::fillPSetDescription(transDesc,
0714                                         false,
0715                                         "KFFitterForRefitInsideOut",
0716                                         "KFSmootherForRefitInsideOut",
0717                                         "PropagatorWithMaterialForMTD",
0718                                         "alongMomentum",
0719                                         true,
0720                                         "WithTrackAngle",
0721                                         "MuonRecHitBuilder",
0722                                         "MTDRecHitBuilder");
0723   desc.add<edm::ParameterSetDescription>("TrackTransformer", transDesc);
0724   desc.add<double>("estimatorMaxChi2", 500.);
0725   desc.add<double>("estimatorMaxNSigma", 10.);
0726   desc.add<double>("btlChi2Cut", 50.);
0727   desc.add<double>("btlTimeChi2Cut", 10.);
0728   desc.add<double>("etlChi2Cut", 50.);
0729   desc.add<double>("etlTimeChi2Cut", 10.);
0730   desc.add<bool>("useVertex", false);
0731   desc.add<bool>("useSimVertex", false);
0732   desc.add<double>("dZCut", 0.1);
0733   desc.add<double>("bsTimeSpread", 0.2);
0734   descriptions.add("trackExtenderWithMTDBase", desc);
0735 }
0736 
0737 template <class TrackCollection>
0738 template <class H, class T>
0739 void TrackExtenderWithMTDT<TrackCollection>::fillValueMap(edm::Event& iEvent,
0740                                                           const H& handle,
0741                                                           const std::vector<T>& vec,
0742                                                           const edm::EDPutToken& token) const {
0743   auto out = std::make_unique<edm::ValueMap<T>>();
0744   typename edm::ValueMap<T>::Filler filler(*out);
0745   filler.insert(handle, vec.begin(), vec.end());
0746   filler.fill();
0747   iEvent.put(token, std::move(out));
0748 }
0749 
0750 template <class TrackCollection>
0751 void TrackExtenderWithMTDT<TrackCollection>::produce(edm::Event& ev, const edm::EventSetup& es) {
0752   //this produces pieces of the track extra
0753   Traj2TrackHits t2t;
0754 
0755   theTransformer->setServices(es);
0756   TrackingRecHitRefProd hitsRefProd = ev.getRefBeforePut<TrackingRecHitCollection>();
0757   reco::TrackExtraRefProd extrasRefProd = ev.getRefBeforePut<reco::TrackExtraCollection>();
0758 
0759   gtg_ = es.getHandle(gtgToken_);
0760 
0761   auto geo = es.getTransientHandle(dlgeoToken_);
0762 
0763   auto magfield = es.getTransientHandle(magfldToken_);
0764 
0765   builder_ = es.getHandle(builderToken_);
0766   hitbuilder_ = es.getHandle(hitbuilderToken_);
0767 
0768   auto propH = es.getTransientHandle(propToken_);
0769   const Propagator* prop = propH.product();
0770 
0771   auto httopo = es.getTransientHandle(ttopoToken_);
0772   const TrackerTopology& ttopo = *httopo;
0773 
0774   auto output = std::make_unique<TrackCollection>();
0775   auto extras = std::make_unique<reco::TrackExtraCollection>();
0776   auto outhits = std::make_unique<edm::OwnVector<TrackingRecHit>>();
0777 
0778   std::vector<float> btlMatchChi2;
0779   std::vector<float> etlMatchChi2;
0780   std::vector<float> btlMatchTimeChi2;
0781   std::vector<float> etlMatchTimeChi2;
0782   std::vector<int> npixBarrel;
0783   std::vector<int> npixEndcap;
0784   std::vector<float> outermostHitPosition;
0785   std::vector<float> pOrigTrkRaw;
0786   std::vector<float> betaOrigTrkRaw;
0787   std::vector<float> t0OrigTrkRaw;
0788   std::vector<float> sigmat0OrigTrkRaw;
0789   std::vector<float> pathLengthsOrigTrkRaw;
0790   std::vector<float> tmtdOrigTrkRaw;
0791   std::vector<float> sigmatmtdOrigTrkRaw;
0792   std::vector<GlobalPoint> tmtdPosOrigTrkRaw;
0793   std::vector<float> tofpiOrigTrkRaw;
0794   std::vector<float> tofkOrigTrkRaw;
0795   std::vector<float> tofpOrigTrkRaw;
0796   std::vector<float> sigmatofpiOrigTrkRaw;
0797   std::vector<float> sigmatofkOrigTrkRaw;
0798   std::vector<float> sigmatofpOrigTrkRaw;
0799   std::vector<int> assocOrigTrkRaw;
0800 
0801   auto const tracksH = ev.getHandle(tracksToken_);
0802 
0803   const auto& trjtrks = ev.get(trajTrackAToken_);
0804 
0805   //MTD hits DetSet
0806   const auto& hits = ev.get(hitsToken_);
0807 
0808   //beam spot
0809   const auto& bs = ev.get(bsToken_);
0810 
0811   const Vertex* pv = nullptr;
0812   if (useVertex_ && !useSimVertex_) {
0813     auto const& vtxs = ev.get(vtxToken_);
0814     if (!vtxs.empty())
0815       pv = &vtxs[0];
0816   }
0817 
0818   std::unique_ptr<math::XYZTLorentzVectorF> genPV(nullptr);
0819   if (useVertex_ && useSimVertex_) {
0820     const auto& genVtxPosition = ev.get(genVtxPositionToken_);
0821     const auto& genVtxTime = ev.get(genVtxTimeToken_);
0822     genPV = std::make_unique<math::XYZTLorentzVectorF>(
0823         genVtxPosition.x(), genVtxPosition.y(), genVtxPosition.z(), genVtxTime);
0824   }
0825 
0826   float vtxTime = 0.f;
0827   if (useVertex_) {
0828     if (useSimVertex_ && genPV) {
0829       vtxTime = genPV->t();
0830     } else if (pv)
0831       vtxTime = pv->t();  //already in ns
0832   }
0833 
0834   std::vector<unsigned> track_indices;
0835   unsigned itrack = 0;
0836 
0837   for (const auto& trjtrk : trjtrks) {
0838     const Trajectory& trajs = *trjtrk.key;
0839     const reco::TrackRef& track = trjtrk.val;
0840 
0841     LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: extrapolating track " << itrack
0842                                      << " p/pT = " << track->p() << " " << track->pt() << " eta = " << track->eta();
0843     LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: sigma_p = "
0844                                      << sqrt(track->covariance()(0, 0)) * track->p2()
0845                                      << " sigma_p/p = " << sqrt(track->covariance()(0, 0)) * track->p() * 100 << " %";
0846 
0847     float trackVtxTime = 0.f;
0848     if (useVertex_) {
0849       float dz;
0850       if (useSimVertex_)
0851         dz = std::abs(track->dz(math::XYZPoint(*genPV)));
0852       else
0853         dz = std::abs(track->dz(pv->position()));
0854 
0855       if (dz < dzCut_)
0856         trackVtxTime = vtxTime;
0857     }
0858 
0859     reco::TransientTrack ttrack(track, magfield.product(), gtg_);
0860     auto thits = theTransformer->getTransientRecHits(ttrack);
0861     TransientTrackingRecHit::ConstRecHitContainer mtdthits;
0862     MTDHitMatchingInfo mBTL, mETL;
0863 
0864     if (trajs.isValid()) {
0865       // get the outermost trajectory point on the track
0866       TrajectoryStateOnSurface tsos = builder_->build(track).outermostMeasurementState();
0867       TrajectoryStateClosestToBeamLine tscbl;
0868       bool tscbl_status = getTrajectoryStateClosestToBeamLine(trajs, bs, prop, tscbl);
0869 
0870       if (tscbl_status) {
0871         float pmag2 = tscbl.trackStateAtPCA().momentum().mag2();
0872         float pathlength0;
0873         TrackSegments trs0;
0874         trackPathLength(trajs, tscbl, prop, pathlength0, trs0);
0875 
0876         const auto& btlhits = tryBTLLayers(tsos,
0877                                            trajs,
0878                                            pmag2,
0879                                            pathlength0,
0880                                            trs0,
0881                                            hits,
0882                                            geo.product(),
0883                                            magfield.product(),
0884                                            prop,
0885                                            bs,
0886                                            trackVtxTime,
0887                                            trackVtxTime != 0.f,
0888                                            mBTL);
0889         mtdthits.insert(mtdthits.end(), btlhits.begin(), btlhits.end());
0890 
0891         // in the future this should include an intermediate refit before propagating to the ETL
0892         // for now it is ok
0893         const auto& etlhits = tryETLLayers(tsos,
0894                                            trajs,
0895                                            pmag2,
0896                                            pathlength0,
0897                                            trs0,
0898                                            hits,
0899                                            geo.product(),
0900                                            magfield.product(),
0901                                            prop,
0902                                            bs,
0903                                            trackVtxTime,
0904                                            trackVtxTime != 0.f,
0905                                            mETL);
0906         mtdthits.insert(mtdthits.end(), etlhits.begin(), etlhits.end());
0907       }
0908 #ifdef EDM_ML_DEBUG
0909       else {
0910         LogTrace("TrackExtenderWithMTD") << "Failing getTrajectoryStateClosestToBeamLine, no search for hits in MTD!";
0911       }
0912 #endif
0913     }
0914 
0915     auto ordering = checkRecHitsOrdering(thits);
0916     if (ordering == RefitDirection::insideOut) {
0917       thits.insert(thits.end(), mtdthits.begin(), mtdthits.end());
0918     } else {
0919       std::reverse(mtdthits.begin(), mtdthits.end());
0920       mtdthits.insert(mtdthits.end(), thits.begin(), thits.end());
0921       thits.swap(mtdthits);
0922     }
0923 
0924     const auto& trajwithmtd =
0925         mtdthits.empty() ? std::vector<Trajectory>(1, trajs) : theTransformer->transform(ttrack, thits);
0926     float pMap = 0.f, betaMap = 0.f, t0Map = 0.f, sigmat0Map = -1.f, pathLengthMap = -1.f, tmtdMap = 0.f,
0927           sigmatmtdMap = -1.f, tofpiMap = 0.f, tofkMap = 0.f, tofpMap = 0.f, sigmatofpiMap = -1.f, sigmatofkMap = -1.f,
0928           sigmatofpMap = -1.f;
0929     GlobalPoint tmtdPosMap{0., 0., 0.};
0930     int iMap = -1;
0931 
0932     for (const auto& trj : trajwithmtd) {
0933       const auto& thetrj = (updateTraj_ ? trj : trajs);
0934       float pathLength = 0.f, tmtd = 0.f, sigmatmtd = -1.f, tofpi = 0.f, tofk = 0.f, tofp = 0.f, sigmatofpi = -1.f,
0935             sigmatofk = -1.f, sigmatofp = -1.f;
0936       GlobalPoint tmtdPos{0., 0., 0.};
0937       LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: refit track " << itrack << " p/pT = " << track->p()
0938                                        << " " << track->pt() << " eta = " << track->eta();
0939       reco::Track result = buildTrack(track,
0940                                       thetrj,
0941                                       trj,
0942                                       bs,
0943                                       magfield.product(),
0944                                       prop,
0945                                       !trajwithmtd.empty() && !mtdthits.empty(),
0946                                       pathLength,
0947                                       tmtd,
0948                                       sigmatmtd,
0949                                       tmtdPos,
0950                                       tofpi,
0951                                       tofk,
0952                                       tofp,
0953                                       sigmatofpi,
0954                                       sigmatofk,
0955                                       sigmatofp);
0956       if (result.ndof() >= 0) {
0957         /// setup the track extras
0958         reco::TrackExtra::TrajParams trajParams;
0959         reco::TrackExtra::Chi2sFive chi2s;
0960         size_t hitsstart = outhits->size();
0961         if (updatePattern_) {
0962           t2t(trj, *outhits, trajParams, chi2s);  // this fills the output hit collection
0963         } else {
0964           t2t(thetrj, *outhits, trajParams, chi2s);
0965         }
0966         size_t hitsend = outhits->size();
0967         extras->push_back(buildTrackExtra(trj));  // always push back the fully built extra, update by setting in track
0968         extras->back().setHits(hitsRefProd, hitsstart, hitsend - hitsstart);
0969         extras->back().setTrajParams(trajParams, chi2s);
0970         //create the track
0971         output->push_back(result);
0972         btlMatchChi2.push_back(mBTL.hit ? mBTL.estChi2 : -1.f);
0973         etlMatchChi2.push_back(mETL.hit ? mETL.estChi2 : -1.f);
0974         btlMatchTimeChi2.push_back(mBTL.hit ? mBTL.timeChi2 : -1.f);
0975         etlMatchTimeChi2.push_back(mETL.hit ? mETL.timeChi2 : -1.f);
0976         pathLengthMap = pathLength;
0977         tmtdMap = tmtd;
0978         sigmatmtdMap = sigmatmtd;
0979         tmtdPosMap = tmtdPos;
0980         auto& backtrack = output->back();
0981         iMap = output->size() - 1;
0982         pMap = backtrack.p();
0983         betaMap = backtrack.beta();
0984         t0Map = backtrack.t0();
0985         sigmat0Map = std::copysign(std::sqrt(std::abs(backtrack.covt0t0())), backtrack.covt0t0());
0986         tofpiMap = tofpi;
0987         tofkMap = tofk;
0988         tofpMap = tofp;
0989         sigmatofpiMap = sigmatofpi;
0990         sigmatofkMap = sigmatofk;
0991         sigmatofpMap = sigmatofp;
0992         reco::TrackExtraRef extraRef(extrasRefProd, extras->size() - 1);
0993         backtrack.setExtra((updateExtra_ ? extraRef : track->extra()));
0994         for (unsigned ihit = hitsstart; ihit < hitsend; ++ihit) {
0995           backtrack.appendHitPattern((*outhits)[ihit], ttopo);
0996         }
0997 #ifdef EDM_ML_DEBUG
0998         LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: hit pattern of refitted track";
0999         for (int i = 0; i < backtrack.hitPattern().numberOfAllHits(reco::HitPattern::TRACK_HITS); i++) {
1000           backtrack.hitPattern().printHitPattern(reco::HitPattern::TRACK_HITS, i, std::cout);
1001         }
1002         LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: missing hit pattern of refitted track";
1003         for (int i = 0; i < backtrack.hitPattern().numberOfAllHits(reco::HitPattern::MISSING_INNER_HITS); i++) {
1004           backtrack.hitPattern().printHitPattern(reco::HitPattern::MISSING_INNER_HITS, i, std::cout);
1005         }
1006 #endif
1007         npixBarrel.push_back(backtrack.hitPattern().numberOfValidPixelBarrelHits());
1008         npixEndcap.push_back(backtrack.hitPattern().numberOfValidPixelEndcapHits());
1009 
1010         if (mBTL.hit || mETL.hit) {
1011           outermostHitPosition.push_back(
1012               mBTL.hit ? (float)(*track).outerRadius()
1013                        : (float)(*track).outerZ());  // save R of the outermost hit for BTL, z for ETL.
1014         } else {
1015           outermostHitPosition.push_back(std::abs(track->eta()) < trackMaxBtlEta_ ? (float)(*track).outerRadius()
1016                                                                                   : (float)(*track).outerZ());
1017         }
1018 
1019         LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: tmtd " << tmtdMap << " +/- " << sigmatmtdMap
1020                                          << " t0 " << t0Map << " +/- " << sigmat0Map << " tof pi/K/p " << tofpiMap
1021                                          << "+/-" << fmt::format("{:0.2g}", sigmatofpiMap) << " ("
1022                                          << fmt::format("{:0.2g}", sigmatofpiMap / tofpiMap * 100) << "%) " << tofkMap
1023                                          << "+/-" << fmt::format("{:0.2g}", sigmatofkMap) << " ("
1024                                          << fmt::format("{:0.2g}", sigmatofkMap / tofkMap * 100) << "%) " << tofpMap
1025                                          << "+/-" << fmt::format("{:0.2g}", sigmatofpMap) << " ("
1026                                          << fmt::format("{:0.2g}", sigmatofpMap / tofpMap * 100) << "%) ";
1027       } else {
1028         LogTrace("TrackExtenderWithMTD") << "Error in the MTD track refitting. This should not happen";
1029       }
1030     }
1031 
1032     pOrigTrkRaw.push_back(pMap);
1033     betaOrigTrkRaw.push_back(betaMap);
1034     t0OrigTrkRaw.push_back(t0Map);
1035     sigmat0OrigTrkRaw.push_back(sigmat0Map);
1036     pathLengthsOrigTrkRaw.push_back(pathLengthMap);
1037     tmtdOrigTrkRaw.push_back(tmtdMap);
1038     sigmatmtdOrigTrkRaw.push_back(sigmatmtdMap);
1039     tmtdPosOrigTrkRaw.push_back(tmtdPosMap);
1040     tofpiOrigTrkRaw.push_back(tofpiMap);
1041     tofkOrigTrkRaw.push_back(tofkMap);
1042     tofpOrigTrkRaw.push_back(tofpMap);
1043     sigmatofpiOrigTrkRaw.push_back(sigmatofpiMap);
1044     sigmatofkOrigTrkRaw.push_back(sigmatofkMap);
1045     sigmatofpOrigTrkRaw.push_back(sigmatofpMap);
1046     assocOrigTrkRaw.push_back(iMap);
1047 
1048     if (iMap == -1) {
1049       btlMatchChi2.push_back(-1.f);
1050       etlMatchChi2.push_back(-1.f);
1051       btlMatchTimeChi2.push_back(-1.f);
1052       etlMatchTimeChi2.push_back(-1.f);
1053       npixBarrel.push_back(-1.f);
1054       npixEndcap.push_back(-1.f);
1055       outermostHitPosition.push_back(0.);
1056     }
1057 
1058     ++itrack;
1059   }
1060 
1061   auto outTrksHandle = ev.put(std::move(output));
1062   ev.put(std::move(extras));
1063   ev.put(std::move(outhits));
1064 
1065   fillValueMap(ev, tracksH, btlMatchChi2, btlMatchChi2Token_);
1066   fillValueMap(ev, tracksH, etlMatchChi2, etlMatchChi2Token_);
1067   fillValueMap(ev, tracksH, btlMatchTimeChi2, btlMatchTimeChi2Token_);
1068   fillValueMap(ev, tracksH, etlMatchTimeChi2, etlMatchTimeChi2Token_);
1069   fillValueMap(ev, tracksH, npixBarrel, npixBarrelToken_);
1070   fillValueMap(ev, tracksH, npixEndcap, npixEndcapToken_);
1071   fillValueMap(ev, tracksH, outermostHitPosition, outermostHitPositionToken_);
1072   fillValueMap(ev, tracksH, pOrigTrkRaw, pOrigTrkToken_);
1073   fillValueMap(ev, tracksH, betaOrigTrkRaw, betaOrigTrkToken_);
1074   fillValueMap(ev, tracksH, t0OrigTrkRaw, t0OrigTrkToken_);
1075   fillValueMap(ev, tracksH, sigmat0OrigTrkRaw, sigmat0OrigTrkToken_);
1076   fillValueMap(ev, tracksH, pathLengthsOrigTrkRaw, pathLengthOrigTrkToken_);
1077   fillValueMap(ev, tracksH, tmtdOrigTrkRaw, tmtdOrigTrkToken_);
1078   fillValueMap(ev, tracksH, sigmatmtdOrigTrkRaw, sigmatmtdOrigTrkToken_);
1079   fillValueMap(ev, tracksH, tmtdPosOrigTrkRaw, tmtdPosOrigTrkToken_);
1080   fillValueMap(ev, tracksH, tofpiOrigTrkRaw, tofpiOrigTrkToken_);
1081   fillValueMap(ev, tracksH, tofkOrigTrkRaw, tofkOrigTrkToken_);
1082   fillValueMap(ev, tracksH, tofpOrigTrkRaw, tofpOrigTrkToken_);
1083   fillValueMap(ev, tracksH, sigmatofpiOrigTrkRaw, sigmatofpiOrigTrkToken_);
1084   fillValueMap(ev, tracksH, sigmatofkOrigTrkRaw, sigmatofkOrigTrkToken_);
1085   fillValueMap(ev, tracksH, sigmatofpOrigTrkRaw, sigmatofpOrigTrkToken_);
1086   fillValueMap(ev, tracksH, assocOrigTrkRaw, assocOrigTrkToken_);
1087 }
1088 
1089 namespace {
1090   bool cmp_for_detset(const unsigned one, const unsigned two) { return one < two; };
1091 
1092   void find_hits_in_dets(const MTDTrackingDetSetVector& hits,
1093                          const Trajectory& traj,
1094                          const DetLayer* layer,
1095                          const TrajectoryStateOnSurface& tsos,
1096                          const float pmag2,
1097                          const float pathlength0,
1098                          const TrackSegments& trs0,
1099                          const float vtxTime,
1100                          const reco::BeamSpot& bs,
1101                          const float bsTimeSpread,
1102                          const Propagator* prop,
1103                          const MeasurementEstimator* estimator,
1104                          bool useVtxConstraint,
1105                          std::set<MTDHitMatchingInfo>& out) {
1106     pair<bool, TrajectoryStateOnSurface> comp = layer->compatible(tsos, *prop, *estimator);
1107     if (comp.first) {
1108       const vector<DetLayer::DetWithState> compDets = layer->compatibleDets(tsos, *prop, *estimator);
1109       LogTrace("TrackExtenderWithMTD") << "Hit search: Compatible dets " << compDets.size();
1110       if (!compDets.empty()) {
1111         for (const auto& detWithState : compDets) {
1112           auto range = hits.equal_range(detWithState.first->geographicalId(), cmp_for_detset);
1113           if (range.first == range.second) {
1114             LogTrace("TrackExtenderWithMTD")
1115                 << "Hit search: no hit in DetId " << detWithState.first->geographicalId().rawId();
1116             continue;
1117           }
1118 
1119           auto pl = prop->propagateWithPath(tsos, detWithState.second.surface());
1120           if (pl.second == 0.) {
1121             LogTrace("TrackExtenderWithMTD")
1122                 << "Hit search: no propagation to DetId " << detWithState.first->geographicalId().rawId();
1123             continue;
1124           }
1125 
1126           const float t_vtx = useVtxConstraint ? vtxTime : 0.f;
1127 
1128           constexpr float vtx_res = 0.008f;
1129           const float t_vtx_err = useVtxConstraint ? vtx_res : bsTimeSpread;
1130 
1131           float lastpmag2 = trs0.getSegmentPathAndMom2(0).second;
1132 
1133           for (auto detitr = range.first; detitr != range.second; ++detitr) {
1134             for (const auto& hit : *detitr) {
1135               auto est = estimator->estimate(detWithState.second, hit);
1136               if (!est.first) {
1137                 LogTrace("TrackExtenderWithMTD")
1138                     << "Hit search: no compatible estimate in DetId " << detWithState.first->geographicalId().rawId()
1139                     << " for hit at pos (" << std::fixed << std::setw(14) << hit.globalPosition().x() << ","
1140                     << std::fixed << std::setw(14) << hit.globalPosition().y() << "," << std::fixed << std::setw(14)
1141                     << hit.globalPosition().z() << ")";
1142                 continue;
1143               }
1144 
1145               LogTrace("TrackExtenderWithMTD")
1146                   << "Hit search: spatial compatibility DetId " << detWithState.first->geographicalId().rawId()
1147                   << " TSOS dx/dy " << std::fixed << std::setw(14)
1148                   << std::sqrt(detWithState.second.localError().positionError().xx()) << " " << std::fixed
1149                   << std::setw(14) << std::sqrt(detWithState.second.localError().positionError().yy()) << " hit dx/dy "
1150                   << std::fixed << std::setw(14) << std::sqrt(hit.localPositionError().xx()) << " " << std::fixed
1151                   << std::setw(14) << std::sqrt(hit.localPositionError().yy()) << " chi2 " << std::fixed
1152                   << std::setw(14) << est.second;
1153 
1154               TrackTofPidInfo tof = computeTrackTofPidInfo(lastpmag2,
1155                                                            std::abs(pl.second),
1156                                                            trs0,
1157                                                            hit.time(),
1158                                                            hit.timeError(),
1159                                                            t_vtx,
1160                                                            t_vtx_err,  //put vtx error by hand for the moment
1161                                                            false,
1162                                                            TofCalc::kMixd,
1163                                                            SigmaTofCalc::kMixd);
1164               MTDHitMatchingInfo mi;
1165               mi.hit = &hit;
1166               mi.estChi2 = est.second;
1167               mi.timeChi2 = tof.dtchi2_best;  //use the chi2 for the best matching hypothesis
1168 
1169               out.insert(mi);
1170             }
1171           }
1172         }
1173       }
1174     }
1175   }
1176 }  // namespace
1177 
1178 template <class TrackCollection>
1179 TransientTrackingRecHit::ConstRecHitContainer TrackExtenderWithMTDT<TrackCollection>::tryBTLLayers(
1180     const TrajectoryStateOnSurface& tsos,
1181     const Trajectory& traj,
1182     const float pmag2,
1183     const float pathlength0,
1184     const TrackSegments& trs0,
1185     const MTDTrackingDetSetVector& hits,
1186     const MTDDetLayerGeometry* geo,
1187     const MagneticField* field,
1188     const Propagator* prop,
1189     const reco::BeamSpot& bs,
1190     const float vtxTime,
1191     const bool matchVertex,
1192     MTDHitMatchingInfo& bestHit) const {
1193   const vector<const DetLayer*>& layers = geo->allBTLLayers();
1194 
1195   TransientTrackingRecHit::ConstRecHitContainer output;
1196   bestHit = MTDHitMatchingInfo();
1197   for (const DetLayer* ilay : layers) {
1198     LogTrace("TrackExtenderWithMTD") << "Hit search: BTL layer at R= "
1199                                      << static_cast<const BarrelDetLayer*>(ilay)->specificSurface().radius();
1200 
1201     fillMatchingHits(ilay, tsos, traj, pmag2, pathlength0, trs0, hits, prop, bs, vtxTime, matchVertex, output, bestHit);
1202   }
1203 
1204   return output;
1205 }
1206 
1207 template <class TrackCollection>
1208 TransientTrackingRecHit::ConstRecHitContainer TrackExtenderWithMTDT<TrackCollection>::tryETLLayers(
1209     const TrajectoryStateOnSurface& tsos,
1210     const Trajectory& traj,
1211     const float pmag2,
1212     const float pathlength0,
1213     const TrackSegments& trs0,
1214     const MTDTrackingDetSetVector& hits,
1215     const MTDDetLayerGeometry* geo,
1216     const MagneticField* field,
1217     const Propagator* prop,
1218     const reco::BeamSpot& bs,
1219     const float vtxTime,
1220     const bool matchVertex,
1221     MTDHitMatchingInfo& bestHit) const {
1222   const vector<const DetLayer*>& layers = geo->allETLLayers();
1223 
1224   TransientTrackingRecHit::ConstRecHitContainer output;
1225   bestHit = MTDHitMatchingInfo();
1226   for (const DetLayer* ilay : layers) {
1227     const BoundDisk& disk = static_cast<const ForwardDetLayer*>(ilay)->specificSurface();
1228     const float diskZ = disk.position().z();
1229 
1230     if (tsos.globalPosition().z() * diskZ < 0)
1231       continue;  // only propagate to the disk that's on the same side
1232 
1233     LogTrace("TrackExtenderWithMTD") << "Hit search: ETL disk at Z = " << diskZ;
1234 
1235     fillMatchingHits(ilay, tsos, traj, pmag2, pathlength0, trs0, hits, prop, bs, vtxTime, matchVertex, output, bestHit);
1236   }
1237 
1238   // the ETL hits order must be from the innermost to the outermost
1239 
1240   if (output.size() == 2) {
1241     if (std::abs(output[0]->globalPosition().z()) > std::abs(output[1]->globalPosition().z())) {
1242       std::reverse(output.begin(), output.end());
1243     }
1244   }
1245   return output;
1246 }
1247 
1248 template <class TrackCollection>
1249 void TrackExtenderWithMTDT<TrackCollection>::fillMatchingHits(const DetLayer* ilay,
1250                                                               const TrajectoryStateOnSurface& tsos,
1251                                                               const Trajectory& traj,
1252                                                               const float pmag2,
1253                                                               const float pathlength0,
1254                                                               const TrackSegments& trs0,
1255                                                               const MTDTrackingDetSetVector& hits,
1256                                                               const Propagator* prop,
1257                                                               const reco::BeamSpot& bs,
1258                                                               const float& vtxTime,
1259                                                               const bool matchVertex,
1260                                                               TransientTrackingRecHit::ConstRecHitContainer& output,
1261                                                               MTDHitMatchingInfo& bestHit) const {
1262   std::set<MTDHitMatchingInfo> hitsInLayer;
1263   bool hitMatched = false;
1264 
1265   using namespace std::placeholders;
1266   auto find_hits = std::bind(find_hits_in_dets,
1267                              std::cref(hits),
1268                              std::cref(traj),
1269                              ilay,
1270                              std::cref(tsos),
1271                              pmag2,
1272                              pathlength0,
1273                              trs0,
1274                              _1,
1275                              std::cref(bs),
1276                              bsTimeSpread_,
1277                              prop,
1278                              theEstimator.get(),
1279                              _2,
1280                              std::ref(hitsInLayer));
1281 
1282   if (useVertex_ && matchVertex)
1283     find_hits(vtxTime, true);
1284   else
1285     find_hits(0, false);
1286 
1287   float spaceChi2Cut = ilay->isBarrel() ? btlChi2Cut_ : etlChi2Cut_;
1288   float timeChi2Cut = ilay->isBarrel() ? btlTimeChi2Cut_ : etlTimeChi2Cut_;
1289 
1290   //just take the first hit because the hits are sorted on their matching quality
1291   if (!hitsInLayer.empty()) {
1292     //check hits to pass minimum quality matching requirements
1293     auto const& firstHit = *hitsInLayer.begin();
1294     LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: matching trial 1: estChi2= " << firstHit.estChi2
1295                                      << " timeChi2= " << firstHit.timeChi2;
1296     if (firstHit.estChi2 < spaceChi2Cut && firstHit.timeChi2 < timeChi2Cut) {
1297       hitMatched = true;
1298       output.push_back(hitbuilder_->build(firstHit.hit));
1299       if (firstHit < bestHit)
1300         bestHit = firstHit;
1301     }
1302   }
1303 
1304   if (useVertex_ && matchVertex && !hitMatched) {
1305     //try a second search with beamspot hypothesis
1306     hitsInLayer.clear();
1307     find_hits(0, false);
1308     if (!hitsInLayer.empty()) {
1309       auto const& firstHit = *hitsInLayer.begin();
1310       LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: matching trial 2: estChi2= " << firstHit.estChi2
1311                                        << " timeChi2= " << firstHit.timeChi2;
1312       if (firstHit.timeChi2 < timeChi2Cut) {
1313         if (firstHit.estChi2 < spaceChi2Cut) {
1314           hitMatched = true;
1315           output.push_back(hitbuilder_->build(firstHit.hit));
1316           if (firstHit < bestHit)
1317             bestHit = firstHit;
1318         }
1319       }
1320     }
1321   }
1322 
1323 #ifdef EDM_ML_DEBUG
1324   if (hitMatched) {
1325     LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: matched hit with time: " << bestHit.hit->time()
1326                                      << " +/- " << bestHit.hit->timeError();
1327   } else {
1328     LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: no matched hit";
1329   }
1330 #endif
1331 }
1332 
1333 //below is unfortunately ripped from other places but
1334 //since track producer doesn't know about MTD we have to do this
1335 template <class TrackCollection>
1336 reco::Track TrackExtenderWithMTDT<TrackCollection>::buildTrack(const reco::TrackRef& orig,
1337                                                                const Trajectory& traj,
1338                                                                const Trajectory& trajWithMtd,
1339                                                                const reco::BeamSpot& bs,
1340                                                                const MagneticField* field,
1341                                                                const Propagator* thePropagator,
1342                                                                bool hasMTD,
1343                                                                float& pathLengthOut,
1344                                                                float& tmtdOut,
1345                                                                float& sigmatmtdOut,
1346                                                                GlobalPoint& tmtdPosOut,
1347                                                                float& tofpi,
1348                                                                float& tofk,
1349                                                                float& tofp,
1350                                                                float& sigmatofpi,
1351                                                                float& sigmatofk,
1352                                                                float& sigmatofp) const {
1353   TrajectoryStateClosestToBeamLine tscbl;
1354   bool tsbcl_status = getTrajectoryStateClosestToBeamLine(traj, bs, thePropagator, tscbl);
1355 
1356   if (!tsbcl_status)
1357     return reco::Track();
1358 
1359   GlobalPoint v = tscbl.trackStateAtPCA().position();
1360   math::XYZPoint pos(v.x(), v.y(), v.z());
1361   GlobalVector p = tscbl.trackStateAtPCA().momentum();
1362   math::XYZVector mom(p.x(), p.y(), p.z());
1363 
1364   int ndof = traj.ndof();
1365 
1366   float t0 = 0.f;
1367   float covt0t0 = -1.f;
1368   pathLengthOut = -1.f;  // if there is no MTD flag the pathlength with -1
1369   tmtdOut = 0.f;
1370   sigmatmtdOut = -1.f;
1371   float betaOut = 0.f;
1372   float covbetabeta = -1.f;
1373 
1374   auto routput = [&]() {
1375     return reco::Track(traj.chiSquared(),
1376                        int(ndof),
1377                        pos,
1378                        mom,
1379                        tscbl.trackStateAtPCA().charge(),
1380                        tscbl.trackStateAtPCA().curvilinearError(),
1381                        orig->algo(),
1382                        reco::TrackBase::undefQuality,
1383                        t0,
1384                        betaOut,
1385                        covt0t0,
1386                        covbetabeta);
1387   };
1388 
1389   //compute path length for time backpropagation, using first MTD hit for the momentum
1390   if (hasMTD) {
1391     float pathlength;
1392     TrackSegments trs;
1393     bool validpropagation = trackPathLength(trajWithMtd, bs, thePropagator, pathlength, trs);
1394     float thit = 0.f;
1395     float thiterror = -1.f;
1396     GlobalPoint thitpos{0., 0., 0.};
1397     bool validmtd = false;
1398 
1399     if (!validpropagation) {
1400       return routput();
1401     }
1402 
1403     uint32_t ihitcount(0), ietlcount(0);
1404     for (auto const& hit : trajWithMtd.measurements()) {
1405       if (hit.recHit()->geographicalId().det() == DetId::Forward &&
1406           ForwardSubdetector(hit.recHit()->geographicalId().subdetId()) == FastTime) {
1407         ihitcount++;
1408         if (MTDDetId(hit.recHit()->geographicalId()).mtdSubDetector() == MTDDetId::MTDType::ETL) {
1409           ietlcount++;
1410         }
1411       }
1412     }
1413 
1414     LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: selected #hits " << ihitcount << " from ETL "
1415                                      << ietlcount;
1416 
1417     auto ihit1 = trajWithMtd.measurements().cbegin();
1418     if (ihitcount == 1) {
1419       const MTDTrackingRecHit* mtdhit = static_cast<const MTDTrackingRecHit*>((*ihit1).recHit()->hit());
1420       thit = mtdhit->time();
1421       thiterror = mtdhit->timeError();
1422       thitpos = mtdhit->globalPosition();
1423       validmtd = true;
1424     } else if (ihitcount == 2 && ietlcount == 2) {
1425       std::pair<float, float> lastStep = trs.getSegmentPathAndMom2(0);
1426       float etlpathlength = std::abs(lastStep.first * c_cm_ns);
1427       //
1428       // The information of the two ETL hits is combined and attributed to the innermost hit
1429       //
1430       if (etlpathlength == 0.f) {
1431         validpropagation = false;
1432       } else {
1433         pathlength -= etlpathlength;
1434         trs.removeFirstSegment();
1435         const MTDTrackingRecHit* mtdhit1 = static_cast<const MTDTrackingRecHit*>((*ihit1).recHit()->hit());
1436         const MTDTrackingRecHit* mtdhit2 = static_cast<const MTDTrackingRecHit*>((*(ihit1 + 1)).recHit()->hit());
1437         TrackTofPidInfo tofInfo = computeTrackTofPidInfo(
1438             lastStep.second, etlpathlength, trs, mtdhit1->time(), mtdhit1->timeError(), 0.f, 0.f, true, TofCalc::kCost);
1439         //
1440         // Protect against incompatible times
1441         //
1442         float err1 = tofInfo.dterror2;
1443         float err2 = mtdhit2->timeError() * mtdhit2->timeError();
1444         if (cms_rounding::roundIfNear0(err1) == 0.f || cms_rounding::roundIfNear0(err2) == 0.f) {
1445           edm::LogError("TrackExtenderWithMTD")
1446               << "MTD tracking hits with zero time uncertainty: " << err1 << " " << err2;
1447         } else {
1448           if ((tofInfo.dt - mtdhit2->time()) * (tofInfo.dt - mtdhit2->time()) < (err1 + err2) * etlTimeChi2Cut_) {
1449             //
1450             // Subtract the ETL time of flight from the outermost measurement, and combine it in a weighted average with the innermost
1451             // the mass ambiguity related uncertainty on the time of flight is added as an additional uncertainty
1452             //
1453             err1 = 1.f / err1;
1454             err2 = 1.f / err2;
1455             thiterror = 1.f / (err1 + err2);
1456             thit = (tofInfo.dt * err1 + mtdhit2->time() * err2) * thiterror;
1457             thiterror = std::sqrt(thiterror);
1458             thitpos = mtdhit2->globalPosition();
1459             LogTrace("TrackExtenderWithMTD")
1460                 << "TrackExtenderWithMTD: p trk = " << p.mag() << " ETL hits times/errors: 1) " << mtdhit1->time()
1461                 << " +/- " << mtdhit1->timeError() << " , 2) " << mtdhit2->time() << " +/- " << mtdhit2->timeError()
1462                 << " extrapolated time1: " << tofInfo.dt << " +/- " << tofInfo.dterror << " average = " << thit
1463                 << " +/- " << thiterror << "\n    hit1 pos: " << mtdhit1->globalPosition()
1464                 << " hit2 pos: " << mtdhit2->globalPosition() << " etl path length " << etlpathlength << std::endl;
1465             validmtd = true;
1466           } else {
1467             // if back extrapolated time of the outermost measurement not compatible with the innermost, keep the one with smallest error
1468             if (err1 <= err2) {
1469               thit = tofInfo.dt;
1470               thiterror = tofInfo.dterror;
1471               validmtd = true;
1472             } else {
1473               thit = mtdhit2->time();
1474               thiterror = mtdhit2->timeError();
1475               validmtd = true;
1476             }
1477           }
1478         }
1479       }
1480     } else {
1481       edm::LogInfo("TrackExtenderWithMTD")
1482           << "MTD hits #" << ihitcount << "ETL hits #" << ietlcount << " anomalous pattern, skipping...";
1483     }
1484 
1485     if (validmtd && validpropagation) {
1486       //here add the PID uncertainty for later use in the 1st step of 4D vtx reconstruction
1487       TrackTofPidInfo tofInfo = computeTrackTofPidInfo(
1488           p.mag2(), pathlength, trs, thit, thiterror, 0.f, 0.f, true, TofCalc::kSegm, SigmaTofCalc::kCost);
1489 
1490       pathLengthOut = pathlength;  // set path length if we've got a timing hit
1491       tmtdOut = thit;
1492       sigmatmtdOut = thiterror;
1493       tmtdPosOut = thitpos;
1494       t0 = tofInfo.dt;
1495       covt0t0 = tofInfo.dterror2;
1496       betaOut = tofInfo.beta_pi;
1497       covbetabeta = tofInfo.betaerror * tofInfo.betaerror;
1498       tofpi = tofInfo.dt_pi;
1499       tofk = tofInfo.dt_k;
1500       tofp = tofInfo.dt_p;
1501       sigmatofpi = tofInfo.sigma_dt_pi;
1502       sigmatofk = tofInfo.sigma_dt_k;
1503       sigmatofp = tofInfo.sigma_dt_p;
1504     }
1505   }
1506 
1507   return routput();
1508 }
1509 
1510 template <class TrackCollection>
1511 reco::TrackExtra TrackExtenderWithMTDT<TrackCollection>::buildTrackExtra(const Trajectory& trajectory) const {
1512   static const string metname = "TrackExtenderWithMTD";
1513 
1514   const Trajectory::RecHitContainer transRecHits = trajectory.recHits();
1515 
1516   // put the collection of TrackingRecHit in the event
1517 
1518   // sets the outermost and innermost TSOSs
1519   // ToDo: validation for track states with MTD
1520   TrajectoryStateOnSurface outerTSOS;
1521   TrajectoryStateOnSurface innerTSOS;
1522   unsigned int innerId = 0, outerId = 0;
1523   TrajectoryMeasurement::ConstRecHitPointer outerRecHit;
1524   DetId outerDetId;
1525 
1526   if (trajectory.direction() == alongMomentum) {
1527     LogTrace(metname) << "alongMomentum";
1528     outerTSOS = trajectory.lastMeasurement().updatedState();
1529     innerTSOS = trajectory.firstMeasurement().updatedState();
1530     outerId = trajectory.lastMeasurement().recHit()->geographicalId().rawId();
1531     innerId = trajectory.firstMeasurement().recHit()->geographicalId().rawId();
1532     outerRecHit = trajectory.lastMeasurement().recHit();
1533     outerDetId = trajectory.lastMeasurement().recHit()->geographicalId();
1534   } else if (trajectory.direction() == oppositeToMomentum) {
1535     LogTrace(metname) << "oppositeToMomentum";
1536     outerTSOS = trajectory.firstMeasurement().updatedState();
1537     innerTSOS = trajectory.lastMeasurement().updatedState();
1538     outerId = trajectory.firstMeasurement().recHit()->geographicalId().rawId();
1539     innerId = trajectory.lastMeasurement().recHit()->geographicalId().rawId();
1540     outerRecHit = trajectory.firstMeasurement().recHit();
1541     outerDetId = trajectory.firstMeasurement().recHit()->geographicalId();
1542   } else
1543     LogError(metname) << "Wrong propagation direction!";
1544 
1545   const GeomDet* outerDet = gtg_->idToDet(outerDetId);
1546   GlobalPoint outerTSOSPos = outerTSOS.globalParameters().position();
1547   bool inside = outerDet->surface().bounds().inside(outerDet->toLocal(outerTSOSPos));
1548 
1549   GlobalPoint hitPos =
1550       (outerRecHit->isValid()) ? outerRecHit->globalPosition() : outerTSOS.globalParameters().position();
1551 
1552   if (!inside) {
1553     LogTrace(metname) << "The Global Muon outerMostMeasurementState is not compatible with the recHit detector!"
1554                       << " Setting outerMost postition to recHit position if recHit isValid: "
1555                       << outerRecHit->isValid();
1556     LogTrace(metname) << "From " << outerTSOSPos << " to " << hitPos;
1557   }
1558 
1559   //build the TrackExtra
1560   GlobalPoint v = (inside) ? outerTSOSPos : hitPos;
1561   GlobalVector p = outerTSOS.globalParameters().momentum();
1562   math::XYZPoint outpos(v.x(), v.y(), v.z());
1563   math::XYZVector outmom(p.x(), p.y(), p.z());
1564 
1565   v = innerTSOS.globalParameters().position();
1566   p = innerTSOS.globalParameters().momentum();
1567   math::XYZPoint inpos(v.x(), v.y(), v.z());
1568   math::XYZVector inmom(p.x(), p.y(), p.z());
1569 
1570   reco::TrackExtra trackExtra(outpos,
1571                               outmom,
1572                               true,
1573                               inpos,
1574                               inmom,
1575                               true,
1576                               outerTSOS.curvilinearError(),
1577                               outerId,
1578                               innerTSOS.curvilinearError(),
1579                               innerId,
1580                               trajectory.direction(),
1581                               trajectory.seedRef());
1582 
1583   return trackExtra;
1584 }
1585 
1586 template <class TrackCollection>
1587 string TrackExtenderWithMTDT<TrackCollection>::dumpLayer(const DetLayer* layer) const {
1588   stringstream output;
1589 
1590   const BoundSurface* sur = nullptr;
1591   const BoundCylinder* bc = nullptr;
1592   const BoundDisk* bd = nullptr;
1593 
1594   sur = &(layer->surface());
1595   if ((bc = dynamic_cast<const BoundCylinder*>(sur))) {
1596     output << "  Cylinder of radius: " << bc->radius() << endl;
1597   } else if ((bd = dynamic_cast<const BoundDisk*>(sur))) {
1598     output << "  Disk at: " << bd->position().z() << endl;
1599   }
1600   return output.str();
1601 }
1602 
1603 //define this as a plug-in
1604 #include <FWCore/Framework/interface/MakerMacros.h>
1605 #include "DataFormats/TrackReco/interface/Track.h"
1606 #include "DataFormats/TrackReco/interface/TrackFwd.h"
1607 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
1608 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
1609 typedef TrackExtenderWithMTDT<reco::TrackCollection> TrackExtenderWithMTD;
1610 
1611 DEFINE_FWK_MODULE(TrackExtenderWithMTD);