56 #include "CLHEP/Units/GlobalPhysicalConstants.h" 70 class MTDHitMatchingInfo {
72 MTDHitMatchingInfo() {
79 inline bool operator<(
const MTDHitMatchingInfo&
m2)
const {
81 constexpr
float chi2_cut = 10.f;
82 constexpr
float low_weight = 3.f;
83 constexpr
float high_weight = 8.f;
84 if (timeChi2 < chi2_cut &&
m2.timeChi2 < chi2_cut)
85 return chi2(low_weight) <
m2.chi2(low_weight);
87 return chi2(high_weight) <
m2.chi2(high_weight);
90 inline float chi2(
float timeWeight = 1.
f)
const {
return estChi2 + timeWeight * timeChi2; }
99 TrackSegments() =
default;
101 inline uint32_t addSegment(
float tPath,
float tMom2) {
102 segmentPathOvc_.emplace_back(tPath *
c_inv);
103 segmentMom2_.emplace_back(tMom2);
106 LogTrace(
"TrackExtenderWithMTD") <<
"addSegment # " << nSegment_ <<
" s = " << tPath
112 inline float computeTof(
float mass_inv2)
const {
114 for (uint32_t iSeg = 0; iSeg < nSegment_; iSeg++) {
115 float gammasq = 1.f + segmentMom2_[iSeg] * mass_inv2;
117 tof += segmentPathOvc_[iSeg] /
beta;
119 LogTrace(
"TrackExtenderWithMTD") <<
" TOF Segment # " << iSeg + 1 <<
" p = " <<
std::sqrt(segmentMom2_[iSeg])
126 inline uint32_t
size()
const {
return nSegment_; }
128 inline uint32_t removeFirstSegment() {
130 segmentPathOvc_.erase(segmentPathOvc_.begin());
131 segmentMom2_.erase(segmentMom2_.begin());
137 inline std::pair<float, float> getSegmentPathAndMom2(uint32_t iSegment)
const {
138 if (iSegment >= nSegment_) {
139 throw cms::Exception(
"TrackExtenderWithMTD") <<
"Requesting non existing track segment #" << iSegment;
141 return std::make_pair(segmentPathOvc_[iSegment], segmentMom2_[iSegment]);
144 uint32_t nSegment_ = 0;
145 std::vector<float> segmentPathOvc_;
146 std::vector<float> segmentMom2_;
149 struct TrackTofPidInfo {
181 enum class TofCalc { kCost = 1, kSegm = 2, kMixd = 3 };
183 const TrackTofPidInfo computeTrackTofPidInfo(
float magp2,
190 bool addPIDError =
true,
191 TofCalc choice = TofCalc::kCost) {
192 constexpr
float m_pi = 0.13957018f;
193 constexpr
float m_pi_inv2 = 1.0f /
m_pi /
m_pi;
194 constexpr
float m_k = 0.493677f;
195 constexpr
float m_k_inv2 = 1.0f / m_k / m_k;
196 constexpr
float m_p = 0.9382720813f;
197 constexpr
float m_p_inv2 = 1.0f / m_p / m_p;
199 TrackTofPidInfo tofpid;
202 tofpid.tmtderror = t_mtderr;
203 tofpid.pathlength = length;
205 auto deltat = [&](
const float mass_inv2,
const float betatmp) {
209 res = tofpid.pathlength / betatmp *
c_inv;
212 res = trs.computeTof(mass_inv2);
215 res = trs.computeTof(mass_inv2) + tofpid.pathlength / betatmp *
c_inv;
221 tofpid.gammasq_pi = 1.f + magp2 * m_pi_inv2;
222 tofpid.beta_pi =
std::sqrt(1.
f - 1.
f / tofpid.gammasq_pi);
223 tofpid.dt_pi = deltat(m_pi_inv2, tofpid.beta_pi);
225 tofpid.gammasq_k = 1.f + magp2 * m_k_inv2;
226 tofpid.beta_k =
std::sqrt(1.
f - 1.
f / tofpid.gammasq_k);
227 tofpid.dt_k = deltat(m_k_inv2, tofpid.beta_k);
229 tofpid.gammasq_p = 1.f + magp2 * m_p_inv2;
230 tofpid.beta_p =
std::sqrt(1.
f - 1.
f / tofpid.gammasq_p);
231 tofpid.dt_p = deltat(m_p_inv2, tofpid.beta_p);
233 tofpid.dt = tofpid.tmtd - tofpid.dt_pi - t_vtx;
234 tofpid.dterror =
sqrt(tofpid.tmtderror * tofpid.tmtderror + t_vtx_err * t_vtx_err);
235 tofpid.betaerror = 0.f;
238 sqrt(tofpid.dterror * tofpid.dterror + (tofpid.dt_p - tofpid.dt_pi) * (tofpid.dt_p - tofpid.dt_pi));
239 tofpid.betaerror = tofpid.beta_p - tofpid.beta_pi;
242 tofpid.dtchi2 = (tofpid.dt * tofpid.dt) / (tofpid.dterror * tofpid.dterror);
244 tofpid.dt_best = tofpid.dt;
245 tofpid.dterror_best = tofpid.dterror;
246 tofpid.dtchi2_best = tofpid.dtchi2;
248 tofpid.prob_pi = -1.f;
249 tofpid.prob_k = -1.f;
250 tofpid.prob_p = -1.f;
254 float chi2_pi = tofpid.dtchi2;
256 (tofpid.tmtd - tofpid.dt_k - t_vtx) * (tofpid.tmtd - tofpid.dt_k - t_vtx) / (tofpid.dterror * tofpid.dterror);
258 (tofpid.tmtd - tofpid.dt_p - t_vtx) * (tofpid.tmtd - tofpid.dt_p - t_vtx) / (tofpid.dterror * tofpid.dterror);
260 float rawprob_pi =
exp(-0.5
f * chi2_pi);
261 float rawprob_k =
exp(-0.5
f * chi2_k);
262 float rawprob_p =
exp(-0.5
f * chi2_p);
263 float normprob = 1.f / (rawprob_pi + rawprob_k + rawprob_p);
265 tofpid.prob_pi = rawprob_pi * normprob;
266 tofpid.prob_k = rawprob_k * normprob;
267 tofpid.prob_p = rawprob_p * normprob;
269 float prob_heavy = 1.f - tofpid.prob_pi;
270 constexpr
float heavy_threshold = 0.75f;
272 if (prob_heavy > heavy_threshold) {
273 if (chi2_k < chi2_p) {
274 tofpid.dt_best = (tofpid.tmtd - tofpid.dt_k - t_vtx);
275 tofpid.dtchi2_best = chi2_k;
277 tofpid.dt_best = (tofpid.tmtd - tofpid.dt_p - t_vtx);
278 tofpid.dtchi2_best = chi2_p;
285 bool getTrajectoryStateClosestToBeamLine(
const Trajectory& traj,
293 if (!stateForProjectionToBeamLineOnSurface.
isValid()) {
294 edm::LogError(
"CannotPropagateToBeamLine") <<
"the state on the closest measurement isnot valid. skipping track.";
301 tscbl = tscblBuilder(stateForProjectionToBeamLine,
bs);
310 TrackSegments& trs) {
313 bool validpropagation =
true;
314 float oldp = traj.
measurements().begin()->updatedState().globalMomentum().mag();
315 float pathlength1 = 0.f;
316 float pathlength2 = 0.f;
320 const auto& propresult = thePropagator->
propagateWithPath(it->updatedState(), (it + 1)->updatedState().surface());
321 float layerpathlength =
std::abs(propresult.second);
322 if (layerpathlength == 0.
f) {
323 validpropagation =
false;
325 pathlength1 += layerpathlength;
326 trs.addSegment(layerpathlength, (it + 1)->updatedState().globalMomentum().
mag2());
328 << std::setw(14) << it->updatedState().globalPosition().perp() <<
" z_i " 329 <<
std::fixed << std::setw(14) << it->updatedState().globalPosition().z()
331 << (it + 1)->updatedState().globalPosition().perp() <<
" z_e " <<
std::fixed 332 << std::setw(14) << (it + 1)->updatedState().globalPosition().z() <<
" p " 333 <<
std::fixed << std::setw(14) << (it + 1)->updatedState().globalMomentum().mag()
335 << (it + 1)->updatedState().globalMomentum().mag() - oldp;
336 oldp = (it + 1)->updatedState().globalMomentum().mag();
344 if (pathlength2 == 0.
f) {
345 validpropagation =
false;
347 pathlength = pathlength1 + pathlength2;
348 trs.addSegment(pathlength2, tscblPCA.momentum().mag2());
350 << std::setw(14) << tscblPCA.position().perp() <<
" z_e " <<
std::fixed 351 << std::setw(14) << tscblPCA.position().z() <<
" p " <<
std::fixed << std::setw(14)
352 << tscblPCA.momentum().mag() <<
" dp " <<
std::fixed << std::setw(14)
353 << tscblPCA.momentum().mag() - oldp;
354 return validpropagation;
361 TrackSegments& trs) {
365 bool tscbl_status = getTrajectoryStateClosestToBeamLine(traj,
bs, thePropagator, tscbl);
370 return trackPathLength(traj, tscbl, thePropagator, pathlength, trs);
375 template <
class TrackCollection>
383 template <
class H,
class T>
394 const TrackSegments&,
401 const bool matchVertex,
402 MTDHitMatchingInfo& bestHit)
const;
408 const TrackSegments&,
415 const bool matchVertex,
416 MTDHitMatchingInfo& bestHit)
const;
418 void fillMatchingHits(
const DetLayer*,
423 const TrackSegments&,
430 MTDHitMatchingInfo&)
const;
439 auto rFirst =
first.mag2();
440 auto rLast =
last.mag2();
446 LogDebug(
"TrackExtenderWithMTD") <<
"Impossible to determine the rechits order" << endl;
465 string dumpLayer(
const DetLayer* layer)
const;
523 template <
class TrackCollection>
529 updateTraj_(iConfig.getParameter<
bool>(
"updateTrackTrajectory")),
530 updateExtra_(iConfig.getParameter<
bool>(
"updateTrackExtra")),
531 updatePattern_(iConfig.getParameter<
bool>(
"updateTrackHitPattern")),
532 mtdRecHitBuilder_(iConfig.getParameter<
std::
string>(
"MTDRecHitBuilder")),
533 propagator_(iConfig.getParameter<
std::
string>(
"Propagator")),
534 transientTrackBuilder_(iConfig.getParameter<
std::
string>(
"TransientTrackBuilder")),
535 estMaxChi2_(iConfig.getParameter<double>(
"estimatorMaxChi2")),
536 estMaxNSigma_(iConfig.getParameter<double>(
"estimatorMaxNSigma")),
537 btlChi2Cut_(iConfig.getParameter<double>(
"btlChi2Cut")),
538 btlTimeChi2Cut_(iConfig.getParameter<double>(
"btlTimeChi2Cut")),
539 etlChi2Cut_(iConfig.getParameter<double>(
"etlChi2Cut")),
540 etlTimeChi2Cut_(iConfig.getParameter<double>(
"etlTimeChi2Cut")),
541 useVertex_(iConfig.getParameter<
bool>(
"useVertex")),
542 useSimVertex_(iConfig.getParameter<
bool>(
"useSimVertex")),
543 dzCut_(iConfig.getParameter<double>(
"dZCut")),
544 bsTimeSpread_(iConfig.getParameter<double>(
"bsTimeSpread")) {
577 gtgToken_ = esConsumes<GlobalTrackingGeometry, GlobalTrackingGeometryRecord>();
578 dlgeoToken_ = esConsumes<MTDDetLayerGeometry, MTDRecoGeometryRecord>();
579 magfldToken_ = esConsumes<MagneticField, IdealMagneticFieldRecord>();
581 ttopoToken_ = esConsumes<TrackerTopology, TrackerTopologyRcd>();
583 produces<edm::OwnVector<TrackingRecHit>>();
584 produces<reco::TrackExtraCollection>();
585 produces<TrackCollection>();
588 template <
class TrackCollection>
598 desc.add<
bool>(
"updateTrackTrajectory",
true);
599 desc.add<
bool>(
"updateTrackExtra",
true);
600 desc.add<
bool>(
"updateTrackHitPattern",
true);
601 desc.add<
std::string>(
"TransientTrackBuilder",
"TransientTrackBuilder");
606 "KFFitterForRefitInsideOut",
607 "KFSmootherForRefitInsideOut",
608 "PropagatorWithMaterialForMTD",
615 desc.add<
double>(
"estimatorMaxChi2", 500.);
616 desc.add<
double>(
"estimatorMaxNSigma", 10.);
617 desc.add<
double>(
"btlChi2Cut", 50.);
618 desc.add<
double>(
"btlTimeChi2Cut", 10.);
619 desc.add<
double>(
"etlChi2Cut", 50.);
620 desc.add<
double>(
"etlTimeChi2Cut", 10.);
621 desc.add<
bool>(
"useVertex",
false);
622 desc.add<
bool>(
"useSimVertex",
false);
623 desc.add<
double>(
"dZCut", 0.1);
624 desc.add<
double>(
"bsTimeSpread", 0.2);
625 descriptions.
add(
"trackExtenderWithMTDBase",
desc);
628 template <
class TrackCollection>
629 template <
class H,
class T>
632 const std::vector<T>& vec,
634 auto out = std::make_unique<edm::ValueMap<T>>();
641 template <
class TrackCollection>
646 theTransformer->setServices(es);
658 hitbuilder_ = es.
getHandle(hitbuilderToken_);
666 auto output = std::make_unique<TrackCollection>();
667 auto extras = std::make_unique<reco::TrackExtraCollection>();
668 auto outhits = std::make_unique<edm::OwnVector<TrackingRecHit>>();
670 std::vector<float> btlMatchChi2;
671 std::vector<float> etlMatchChi2;
672 std::vector<float> btlMatchTimeChi2;
673 std::vector<float> etlMatchTimeChi2;
674 std::vector<int> npixBarrel;
675 std::vector<int> npixEndcap;
676 std::vector<float> pOrigTrkRaw;
677 std::vector<float> betaOrigTrkRaw;
678 std::vector<float> t0OrigTrkRaw;
679 std::vector<float> sigmat0OrigTrkRaw;
680 std::vector<float> pathLengthsOrigTrkRaw;
681 std::vector<float> tmtdOrigTrkRaw;
682 std::vector<float> sigmatmtdOrigTrkRaw;
683 std::vector<float> tofpiOrigTrkRaw;
684 std::vector<float> tofkOrigTrkRaw;
685 std::vector<float> tofpOrigTrkRaw;
686 std::vector<int> assocOrigTrkRaw;
688 auto const tracksH =
ev.getHandle(tracksToken_);
690 const auto& trjtrks =
ev.get(trajTrackAToken_);
693 const auto&
hits =
ev.get(hitsToken_);
696 const auto&
bs =
ev.get(bsToken_);
699 if (useVertex_ && !useSimVertex_) {
700 auto const& vtxs =
ev.get(vtxToken_);
705 std::unique_ptr<math::XYZTLorentzVectorF> genPV(
nullptr);
706 if (useVertex_ && useSimVertex_) {
707 const auto& genVtxPosition =
ev.get(genVtxPositionToken_);
708 const auto& genVtxTime =
ev.get(genVtxTimeToken_);
709 genPV = std::make_unique<math::XYZTLorentzVectorF>(
710 genVtxPosition.x(), genVtxPosition.y(), genVtxPosition.z(), genVtxTime);
715 if (useSimVertex_ && genPV) {
716 vtxTime = genPV->t();
721 std::vector<unsigned> track_indices;
724 for (
const auto& trjtrk : trjtrks) {
728 LogTrace(
"TrackExtenderWithMTD") <<
"TrackExtenderWithMTD: extrapolating track " << itrack
729 <<
" p/pT = " <<
track->p() <<
" " <<
track->pt() <<
" eta = " <<
track->eta();
731 float trackVtxTime = 0.f;
740 trackVtxTime = vtxTime;
744 auto thits = theTransformer->getTransientRecHits(ttrack);
746 MTDHitMatchingInfo mBTL, mETL;
752 bool tscbl_status = getTrajectoryStateClosestToBeamLine(trajs,
bs, prop, tscbl);
758 trackPathLength(trajs, tscbl, prop, pathlength0, trs0);
760 const auto& btlhits = tryBTLLayers(tsos,
773 mtdthits.insert(mtdthits.end(), btlhits.begin(), btlhits.end());
777 const auto& etlhits = tryETLLayers(tsos,
790 mtdthits.insert(mtdthits.end(), etlhits.begin(), etlhits.end());
794 auto ordering = checkRecHitsOrdering(thits);
796 thits.insert(thits.end(), mtdthits.begin(), mtdthits.end());
799 mtdthits.insert(mtdthits.end(), thits.begin(), thits.end());
800 thits.swap(mtdthits);
803 const auto& trajwithmtd =
804 mtdthits.empty() ? std::vector<Trajectory>(1, trajs) : theTransformer->transform(ttrack, thits);
805 float pMap = 0.f, betaMap = 0.f, t0Map = 0.f, sigmat0Map = -1.f, pathLengthMap = -1.f, tmtdMap = 0.f,
806 sigmatmtdMap = -1.f, tofpiMap = 0.f, tofkMap = 0.f, tofpMap = 0.f;
809 for (
const auto& trj : trajwithmtd) {
810 const auto& thetrj = (updateTraj_ ? trj : trajs);
811 float pathLength = 0.f, tmtd = 0.f, sigmatmtd = -1.f, tofpi = 0.f, tofk = 0.f, tofp = 0.f;
812 LogTrace(
"TrackExtenderWithMTD") <<
"TrackExtenderWithMTD: refit track " << itrack <<
" p/pT = " <<
track->p()
813 <<
" " <<
track->pt() <<
" eta = " <<
track->eta();
820 !trajwithmtd.empty() && !mtdthits.empty(),
831 size_t hitsstart = outhits->size();
832 if (updatePattern_) {
833 t2t(trj, *outhits, trajParams, chi2s);
835 t2t(thetrj, *outhits, trajParams, chi2s);
837 size_t hitsend = outhits->size();
838 extras->push_back(buildTrackExtra(trj));
839 extras->back().setHits(hitsRefProd, hitsstart, hitsend - hitsstart);
840 extras->back().setTrajParams(trajParams, chi2s);
843 btlMatchChi2.push_back(mBTL.hit ? mBTL.estChi2 : -1.f);
844 etlMatchChi2.push_back(mETL.hit ? mETL.estChi2 : -1.f);
845 btlMatchTimeChi2.push_back(mBTL.hit ? mBTL.timeChi2 : -1.f);
846 etlMatchTimeChi2.push_back(mETL.hit ? mETL.timeChi2 : -1.f);
847 pathLengthMap = pathLength;
849 sigmatmtdMap = sigmatmtd;
850 auto& backtrack =
output->back();
851 iMap =
output->size() - 1;
852 pMap = backtrack.p();
853 betaMap = backtrack.beta();
854 t0Map = backtrack.t0();
855 sigmat0Map = std::copysign(
std::sqrt(
std::abs(backtrack.covt0t0())), backtrack.covt0t0());
860 backtrack.setExtra((updateExtra_ ? extraRef :
track->extra()));
861 for (
unsigned ihit = hitsstart; ihit < hitsend; ++ihit) {
862 backtrack.appendHitPattern((*outhits)[ihit], ttopo);
864 npixBarrel.push_back(backtrack.hitPattern().numberOfValidPixelBarrelHits());
865 npixEndcap.push_back(backtrack.hitPattern().numberOfValidPixelEndcapHits());
866 LogTrace(
"TrackExtenderWithMTD") <<
"TrackExtenderWithMTD: tmtd " << tmtdMap <<
" +/- " << sigmatmtdMap
867 <<
" t0 " << t0Map <<
" +/- " << sigmat0Map <<
" tof pi/K/p " << tofpiMap
868 <<
" " << tofkMap <<
" " << tofpMap;
870 LogTrace(
"TrackExtenderWithMTD") <<
"Error in the MTD track refitting. This should not happen";
874 pOrigTrkRaw.push_back(pMap);
875 betaOrigTrkRaw.push_back(betaMap);
876 t0OrigTrkRaw.push_back(t0Map);
877 sigmat0OrigTrkRaw.push_back(sigmat0Map);
878 pathLengthsOrigTrkRaw.push_back(pathLengthMap);
879 tmtdOrigTrkRaw.push_back(tmtdMap);
880 sigmatmtdOrigTrkRaw.push_back(sigmatmtdMap);
881 tofpiOrigTrkRaw.push_back(tofpiMap);
882 tofkOrigTrkRaw.push_back(tofkMap);
883 tofpOrigTrkRaw.push_back(tofpMap);
884 assocOrigTrkRaw.push_back(iMap);
887 btlMatchChi2.push_back(-1.
f);
888 etlMatchChi2.push_back(-1.
f);
889 btlMatchTimeChi2.push_back(-1.
f);
890 etlMatchTimeChi2.push_back(-1.
f);
891 npixBarrel.push_back(-1.
f);
892 npixEndcap.push_back(-1.
f);
902 fillValueMap(
ev, tracksH, btlMatchChi2, btlMatchChi2Token_);
903 fillValueMap(
ev, tracksH, etlMatchChi2, etlMatchChi2Token_);
904 fillValueMap(
ev, tracksH, btlMatchTimeChi2, btlMatchTimeChi2Token_);
905 fillValueMap(
ev, tracksH, etlMatchTimeChi2, etlMatchTimeChi2Token_);
906 fillValueMap(
ev, tracksH, npixBarrel, npixBarrelToken_);
907 fillValueMap(
ev, tracksH, npixEndcap, npixEndcapToken_);
908 fillValueMap(
ev, tracksH, pOrigTrkRaw, pOrigTrkToken_);
909 fillValueMap(
ev, tracksH, betaOrigTrkRaw, betaOrigTrkToken_);
910 fillValueMap(
ev, tracksH, t0OrigTrkRaw, t0OrigTrkToken_);
911 fillValueMap(
ev, tracksH, sigmat0OrigTrkRaw, sigmat0OrigTrkToken_);
912 fillValueMap(
ev, tracksH, pathLengthsOrigTrkRaw, pathLengthOrigTrkToken_);
913 fillValueMap(
ev, tracksH, tmtdOrigTrkRaw, tmtdOrigTrkToken_);
914 fillValueMap(
ev, tracksH, sigmatmtdOrigTrkRaw, sigmatmtdOrigTrkToken_);
915 fillValueMap(
ev, tracksH, tofpiOrigTrkRaw, tofpiOrigTrkToken_);
916 fillValueMap(
ev, tracksH, tofkOrigTrkRaw, tofkOrigTrkToken_);
917 fillValueMap(
ev, tracksH, tofpOrigTrkRaw, tofpOrigTrkToken_);
918 fillValueMap(
ev, tracksH, assocOrigTrkRaw, assocOrigTrkToken_);
922 bool cmp_for_detset(
const unsigned one,
const unsigned two) {
return one <
two; };
929 const float pathlength0,
930 const TrackSegments& trs0,
933 const float bsTimeSpread,
936 bool useVtxConstraint,
937 std::set<MTDHitMatchingInfo>&
out) {
938 pair<bool, TrajectoryStateOnSurface>
comp =
layer->compatible(tsos, *prop, *
estimator);
940 const vector<DetLayer::DetWithState> compDets =
layer->compatibleDets(tsos, *prop, *
estimator);
941 LogTrace(
"TrackExtenderWithMTD") <<
"Hit search: Compatible dets " << compDets.size();
942 if (!compDets.empty()) {
943 for (
const auto& detWithState : compDets) {
944 auto range =
hits.equal_range(detWithState.first->geographicalId(), cmp_for_detset);
947 <<
"Hit search: no hit in DetId " << detWithState.first->geographicalId().rawId();
952 if (pl.second == 0.) {
954 <<
"Hit search: no propagation to DetId " << detWithState.first->geographicalId().rawId();
958 const float t_vtx = useVtxConstraint ? vtxTime : 0.f;
960 constexpr
float vtx_res = 0.008f;
961 const float t_vtx_err = useVtxConstraint ? vtx_res : bsTimeSpread;
963 float lastpmag2 = trs0.getSegmentPathAndMom2(0).second;
965 for (
auto detitr =
range.first; detitr !=
range.second; ++detitr) {
966 for (
const auto&
hit : *detitr) {
967 auto est =
estimator->estimate(detWithState.second,
hit);
970 <<
"Hit search: no compatible estimate in DetId " << detWithState.first->geographicalId().rawId()
971 <<
" for hit at pos (" <<
std::fixed << std::setw(14) <<
hit.globalPosition().
x() <<
"," 973 <<
hit.globalPosition().
z() <<
")";
978 <<
"Hit search: spatial compatibility DetId " << detWithState.first->geographicalId().rawId()
979 <<
" TSOS dx/dy " <<
std::fixed << std::setw(14)
981 << std::setw(14) <<
std::sqrt(detWithState.second.localError().positionError().yy()) <<
" hit dx/dy " 984 << std::setw(14) << est.second;
986 TrackTofPidInfo tof = computeTrackTofPidInfo(lastpmag2,
995 MTDHitMatchingInfo mi;
997 mi.estChi2 = est.second;
998 mi.timeChi2 = tof.dtchi2_best;
1009 template <
class TrackCollection>
1014 const float pathlength0,
1015 const TrackSegments& trs0,
1021 const float vtxTime,
1022 const bool matchVertex,
1023 MTDHitMatchingInfo& bestHit)
const {
1027 bestHit = MTDHitMatchingInfo();
1029 fillMatchingHits(ilay, tsos, traj, pmag2, pathlength0, trs0,
hits, prop,
bs, vtxTime, matchVertex,
output, bestHit);
1033 template <
class TrackCollection>
1038 const float pathlength0,
1039 const TrackSegments& trs0,
1045 const float vtxTime,
1046 const bool matchVertex,
1047 MTDHitMatchingInfo& bestHit)
const {
1051 bestHit = MTDHitMatchingInfo();
1054 const float diskZ = disk.position().z();
1059 fillMatchingHits(ilay, tsos, traj, pmag2, pathlength0, trs0,
hits, prop,
bs, vtxTime, matchVertex,
output, bestHit);
1064 if (
output.size() == 2) {
1072 template <
class TrackCollection>
1077 const float pathlength0,
1078 const TrackSegments& trs0,
1082 const float& vtxTime,
1083 const bool matchVertex,
1085 MTDHitMatchingInfo& bestHit)
const {
1086 std::set<MTDHitMatchingInfo> hitsInLayer;
1087 bool hitMatched =
false;
1090 auto find_hits = std::bind(find_hits_in_dets,
1104 std::ref(hitsInLayer));
1106 if (useVertex_ && matchVertex)
1107 find_hits(vtxTime,
true);
1109 find_hits(0,
false);
1111 float spaceChi2Cut = ilay->
isBarrel() ? btlChi2Cut_ : etlChi2Cut_;
1112 float timeChi2Cut = ilay->
isBarrel() ? btlTimeChi2Cut_ : etlTimeChi2Cut_;
1115 if (!hitsInLayer.empty()) {
1117 auto const& firstHit = *hitsInLayer.begin();
1118 LogTrace(
"TrackExtenderWithMTD") <<
"TrackExtenderWithMTD: estChi2= " << firstHit.estChi2
1119 <<
" timeChi2= " << firstHit.timeChi2;
1120 if (firstHit.estChi2 < spaceChi2Cut && firstHit.timeChi2 < timeChi2Cut) {
1122 output.push_back(hitbuilder_->build(firstHit.hit));
1123 if (firstHit < bestHit)
1128 if (useVertex_ && matchVertex && !hitMatched) {
1130 hitsInLayer.clear();
1131 find_hits(0,
false);
1132 if (!hitsInLayer.empty()) {
1133 auto const& firstHit = *hitsInLayer.begin();
1134 LogTrace(
"TrackExtenderWithMTD") <<
"TrackExtenderWithMTD: estChi2= " << firstHit.estChi2
1135 <<
" timeChi2= " << firstHit.timeChi2;
1136 if (firstHit.timeChi2 < timeChi2Cut) {
1137 if (firstHit.estChi2 < spaceChi2Cut) {
1139 output.push_back(hitbuilder_->build(firstHit.hit));
1140 if (firstHit < bestHit)
1150 template <
class TrackCollection>
1158 float& pathLengthOut,
1160 float& sigmatmtdOut,
1163 float& tofp)
const {
1165 bool tsbcl_status = getTrajectoryStateClosestToBeamLine(traj,
bs, thePropagator, tscbl);
1178 float covt0t0 = -1.f;
1179 pathLengthOut = -1.f;
1181 sigmatmtdOut = -1.f;
1182 float betaOut = 0.f;
1183 float covbetabeta = -1.f;
1185 auto routput = [&]() {
1204 bool validpropagation = trackPathLength(trajWithMtd,
bs, thePropagator, pathlength, trs);
1206 float thiterror = -1.f;
1207 bool validmtd =
false;
1209 if (!validpropagation) {
1213 uint32_t ihitcount(0), ietlcount(0);
1218 if (
MTDDetId(
hit.recHit()->geographicalId()).mtdSubDetector() == MTDDetId::MTDType::ETL) {
1225 if (ihitcount == 1) {
1227 thit = mtdhit->
time();
1230 }
else if (ihitcount == 2 && ietlcount == 2) {
1231 std::pair<float, float> lastStep = trs.getSegmentPathAndMom2(0);
1236 if (etlpathlength == 0.
f) {
1237 validpropagation =
false;
1239 pathlength -= etlpathlength;
1240 trs.removeFirstSegment();
1243 TrackTofPidInfo tofInfo = computeTrackTofPidInfo(
1244 lastStep.second, etlpathlength, trs, mtdhit1->
time(), mtdhit1->
timeError(), 0.f, 0.f,
true, TofCalc::kCost);
1248 float err1 = tofInfo.dterror * tofInfo.dterror;
1252 <<
"MTD tracking hits with zero time uncertainty: " << err1 <<
" " << err2;
1254 if ((tofInfo.dt - mtdhit2->
time()) * (tofInfo.dt - mtdhit2->
time()) < (err1 + err2) * etlTimeChi2Cut_) {
1261 thiterror = 1.f / (err1 + err2);
1262 thit = (tofInfo.dt * err1 + mtdhit2->
time() * err2) * thiterror;
1265 <<
"TrackExtenderWithMTD: p trk = " <<
p.mag() <<
" ETL hits times/errors: " << mtdhit1->
time()
1267 <<
" extrapolated time1: " << tofInfo.dt <<
" +/- " << tofInfo.dterror <<
" average = " << thit
1268 <<
" +/- " << thiterror;
1275 <<
"MTD hits #" << ihitcount <<
"ETL hits #" << ietlcount <<
" anomalous pattern, skipping...";
1278 if (validmtd && validpropagation) {
1280 TrackTofPidInfo tofInfo =
1281 computeTrackTofPidInfo(
p.mag2(), pathlength, trs, thit, thiterror, 0.f, 0.f,
true, TofCalc::kSegm);
1282 pathLengthOut = pathlength;
1284 sigmatmtdOut = thiterror;
1286 covt0t0 = tofInfo.dterror * tofInfo.dterror;
1287 betaOut = tofInfo.beta_pi;
1288 covbetabeta = tofInfo.betaerror * tofInfo.betaerror;
1289 tofpi = tofInfo.dt_pi;
1290 tofk = tofInfo.dt_k;
1291 tofp = tofInfo.dt_p;
1298 template <
class TrackCollection>
1300 static const string metname =
"TrackExtenderWithMTD";
1310 unsigned int innerId = 0, outerId = 0;
1333 const GeomDet* outerDet = gtg_->idToDet(outerDetId);
1341 LogTrace(
metname) <<
"The Global Muon outerMostMeasurementState is not compatible with the recHit detector!" 1342 <<
" Setting outerMost postition to recHit position if recHit isValid: " 1343 << outerRecHit->isValid();
1374 template <
class TrackCollection>
1382 sur = &(
layer->surface());
1383 if ((bc = dynamic_cast<const BoundCylinder*>(sur))) {
1384 output <<
" Cylinder of radius: " << bc->radius() << endl;
1385 }
else if ((bd = dynamic_cast<const BoundDisk*>(sur))) {
1386 output <<
" Disk at: " << bd->position().z() << endl;
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
TrackExtenderWithMTDT< reco::TrackCollection > TrackExtenderWithMTD
edm::EDPutToken btlMatchTimeChi2Token_
T getParameter(std::string const &) const
edm::EDPutToken pOrigTrkToken_
edm::ESGetToken< TransientTrackingRecHitBuilder, TransientRecHitRecord > hitbuilderToken_
edm::EDPutToken t0OrigTrkToken_
const CurvilinearTrajectoryError & curvilinearError() const
const float bsTimeSpread_
edm::EDPutToken btlMatchChi2Token_
virtual bool inside(const Local3DPoint &) const =0
Determine if the point is inside the bounds.
edm::EDPutToken tofpOrigTrkToken_
const float btlTimeChi2Cut_
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
GlobalPoint position() const
edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > gtgToken_
const std::string metname
edm::EDGetTokenT< reco::BeamSpot > bsToken_
std::unique_ptr< MeasurementEstimator > theEstimator
TrackExtenderWithMTDT(const ParameterSet &pset)
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > ttopoToken_
const GlobalTrajectoryParameters & globalParameters() const
ParameterSet const & getParameterSet(std::string const &) const
Global3DPoint GlobalPoint
edm::EDPutToken etlMatchChi2Token_
edm::View< TrackType > InputCollection
edm::EDGetTokenT< GlobalPoint > genVtxPositionToken_
edm::RefToBase< TrajectorySeed > seedRef(void) const
Log< level::Error, false > LogError
const SurfaceType & surface() const
edm::EDPutToken tofkOrigTrkToken_
edm::EDPutToken etlMatchTimeChi2Token_
edm::EDPutToken sigmat0OrigTrkToken_
void fillValueMap(edm::Event &iEvent, const H &handle, const std::vector< T > &vec, const edm::EDPutToken &token) const
edm::EDPutToken tmtdOrigTrkToken_
TrajectoryMeasurement const & lastMeasurement() const
edm::EDGetTokenT< MTDTrackingDetSetVector > hitsToken_
Detector identifier base class for the MIP Timing Layer.
TransientTrackingRecHit::ConstRecHitContainer tryETLLayers(const TrajectoryStateOnSurface &, const Trajectory &traj, const float, const float, const TrackSegments &, const MTDTrackingDetSetVector &, const MTDDetLayerGeometry *, const MagneticField *field, const Propagator *prop, const reco::BeamSpot &bs, const float vtxTime, const bool matchVertex, MTDHitMatchingInfo &bestHit) const
DataContainer const & measurements() const
GlobalPoint position() const
string dumpLayer(const DetLayer *layer) const
GlobalVector momentum() const
int ndof(bool bon=true) const
Container::value_type value_type
GlobalPoint globalPosition() const
const float estMaxNSigma_
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > builderToken_
const std::string mtdRecHitBuilder_
void fillMatchingHits(const DetLayer *, const TrajectoryStateOnSurface &, const Trajectory &, const float, const float, const TrackSegments &, const MTDTrackingDetSetVector &, const Propagator *, const reco::BeamSpot &, const float &, const bool, TransientTrackingRecHit::ConstRecHitContainer &, MTDHitMatchingInfo &) const
PropagationDirection const & direction() const
TrackCharge charge() const
ConstRecHitContainer recHits() const
GlobalVector momentum() const
Abs< T >::type abs(const T &t)
#define DEFINE_FWK_MODULE(type)
FTS const & trackStateAtPCA() const
edm::EDPutToken assocOrigTrkToken_
const std::string propagator_
edm::ESGetToken< Propagator, TrackingComponentsRecord > propToken_
virtual std::pair< TrajectoryStateOnSurface, double > propagateWithPath(const FreeTrajectoryState &, const Surface &) const final
constexpr valType roundIfNear0(valType value, double tolerance=1.e-7)
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
edm::ESHandle< TransientTrackingRecHitBuilder > hitbuilder_
const std::vector< const DetLayer * > & allBTLLayers() const
return the BTL DetLayers (barrel), inside-out
std::vector< ConstRecHitPointer > ConstRecHitContainer
const std::string transientTrackBuilder_
Log< level::Info, false > LogInfo
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
ConstRecHitContainer RecHitContainer
constexpr NumType convertMmToCm(NumType millimeters)
const Plane & surface() const
The nominal surface of the GeomDet.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const CurvilinearTrajectoryError & curvilinearError() const
XYZVectorD XYZVector
spatial vector with cartesian internal representation
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
XYZPointD XYZPoint
point in space with cartesian internal representation
reco::TrackExtra buildTrackExtra(const Trajectory &trajectory) const
A 2D TrackerRecHit with time and time error information.
void add(std::string const &label, ParameterSetDescription const &psetDescription)
TrackingRecHit::ConstRecHitPointer ConstRecHitPointer
edm::ESHandle< TransientTrackBuilder > builder_
std::unique_ptr< TrackTransformer > theTransformer
TrajectoryStateOnSurface const & updatedState() const
edm::EDPutToken npixBarrelToken_
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magfldToken_
bool operator<(DTCELinkId const &lhs, DTCELinkId const &rhs)
edm::EDPutToken npixEndcapToken_
TrajectoryMeasurement const & firstMeasurement() const
static int position[264][3]
TrajectoryMeasurement const & closestMeasurement(GlobalPoint) const
edm::EDPutToken betaOrigTrkToken_
FreeTrajectoryState const * freeState(bool withErrors=true) const
ESTransientHandle< T > getTransientHandle(const ESGetToken< T, R > &iToken) const
edm::EDGetTokenT< TrajTrackAssociationCollection > trajTrackAToken_
edm::EDPutToken pathLengthOrigTrkToken_
edm::ESHandle< GlobalTrackingGeometry > gtg_
edm::EDPutToken sigmatmtdOrigTrkToken_
TransientTrackingRecHit::ConstRecHitContainer tryBTLLayers(const TrajectoryStateOnSurface &, const Trajectory &traj, const float, const float, const TrackSegments &, const MTDTrackingDetSetVector &, const MTDDetLayerGeometry *, const MagneticField *field, const Propagator *prop, const reco::BeamSpot &bs, const float vtxTime, const bool matchVertex, MTDHitMatchingInfo &bestHit) const
edm::EDPutToken tofpiOrigTrkToken_
edm::EDGetTokenT< float > genVtxTimeToken_
edm::EDGetTokenT< VertexCollection > vtxToken_
void produce(edm::Event &ev, const edm::EventSetup &es) final
TrackCollection::value_type TrackType
edm::EDGetTokenT< InputCollection > tracksToken_
const std::vector< const DetLayer * > & allETLLayers() const
return the ETL DetLayers (endcap), -Z to +Z
const float etlTimeChi2Cut_
RefitDirection::GeometricalDirection checkRecHitsOrdering(TransientTrackingRecHit::ConstRecHitContainer const &recHits) const
ConstRecHitPointer const & recHit() const
reco::Track buildTrack(const reco::TrackRef &, const Trajectory &, const Trajectory &, const reco::BeamSpot &, const MagneticField *field, const Propagator *prop, bool hasMTD, float &pathLength, float &tmtdOut, float &sigmatmtdOut, float &tofpi, float &tofk, float &tofp) const
edm::ESGetToken< MTDDetLayerGeometry, MTDRecoGeometryRecord > dlgeoToken_
const Bounds & bounds() const