Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-12 23:04:07

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