CMS 3D CMS Logo

TrackExtenderWithMTD.cc
Go to the documentation of this file.
8 
11 
14 
16 
18 
21 
27 
30 
34 
36 
39 
42 
44 
47 
48 #include <sstream>
49 
51 
54 
57 #include "CLHEP/Units/GlobalPhysicalConstants.h"
59 
62 
63 using namespace std;
64 using namespace edm;
65 using namespace reco;
66 
67 namespace {
68  constexpr float c_cm_ns = geant_units::operators::convertMmToCm(CLHEP::c_light); // [mm/ns] -> [cm/ns]
69  constexpr float c_inv = 1.0f / c_cm_ns;
70 
71  class MTDHitMatchingInfo {
72  public:
73  MTDHitMatchingInfo() {
74  hit = nullptr;
77  }
78 
79  //Operator used to sort the hits while performing the matching step at the MTD
80  inline bool operator<(const MTDHitMatchingInfo& m2) const {
81  //only for good matching in time use estChi2, otherwise use mostly time compatibility
82  constexpr float chi2_cut = 10.f;
83  constexpr float low_weight = 3.f;
84  constexpr float high_weight = 8.f;
85  if (timeChi2 < chi2_cut && m2.timeChi2 < chi2_cut)
86  return chi2(low_weight) < m2.chi2(low_weight);
87  else
88  return chi2(high_weight) < m2.chi2(high_weight);
89  }
90 
91  inline float chi2(float timeWeight = 1.f) const { return estChi2 + timeWeight * timeChi2; }
92 
93  const MTDTrackingRecHit* hit;
94  float estChi2;
95  float timeChi2;
96  };
97 
98  class TrackSegments {
99  public:
100  TrackSegments() {
101  sigmaTofs_.reserve(30); // observed upper limit on nSegments
102  };
103 
104  inline uint32_t addSegment(float tPath, float tMom2, float sigmaMom) {
105  segmentPathOvc_.emplace_back(tPath * c_inv);
106  segmentMom2_.emplace_back(tMom2);
107  segmentSigmaMom_.emplace_back(sigmaMom);
108  nSegment_++;
109 
110  LogTrace("TrackExtenderWithMTD") << "addSegment # " << nSegment_ << " s = " << tPath
111  << " p = " << std::sqrt(tMom2) << " sigma_p = " << sigmaMom
112  << " sigma_p/p = " << sigmaMom / std::sqrt(tMom2) * 100 << " %";
113 
114  return nSegment_;
115  }
116 
117  inline float computeTof(float mass_inv2) const {
118  float tof(0.f);
119  for (uint32_t iSeg = 0; iSeg < nSegment_; iSeg++) {
120  float gammasq = 1.f + segmentMom2_[iSeg] * mass_inv2;
121  float beta = std::sqrt(1.f - 1.f / gammasq);
122  tof += segmentPathOvc_[iSeg] / beta;
123 
124  LogTrace("TrackExtenderWithMTD") << " TOF Segment # " << iSeg + 1 << " p = " << std::sqrt(segmentMom2_[iSeg])
125  << " tof = " << tof;
126 
127 #ifdef EDM_ML_DEBUG
128  float sigma_tof = segmentPathOvc_[iSeg] * segmentSigmaMom_[iSeg] /
129  (segmentMom2_[iSeg] * sqrt(segmentMom2_[iSeg] + 1 / mass_inv2) * mass_inv2);
130 
131  LogTrace("TrackExtenderWithMTD") << "TOF Segment # " << iSeg + 1 << std::fixed << std::setw(6)
132  << " tof segment = " << segmentPathOvc_[iSeg] / beta << std::scientific
133  << "+/- " << sigma_tof << std::fixed
134  << "(rel. err. = " << sigma_tof / (segmentPathOvc_[iSeg] / beta) * 100
135  << " %)";
136 #endif
137  }
138 
139  return tof;
140  }
141 
142  inline float computeSigmaTof(float mass_inv2) {
143  float sigmatof = 0.;
144 
145  // remove previously calculated sigmaTofs
146  sigmaTofs_.clear();
147 
148  // compute sigma(tof) on each segment first by propagating sigma(p)
149  // also add diagonal terms to sigmatof
150  float sigma = 0.;
151  for (uint32_t iSeg = 0; iSeg < nSegment_; iSeg++) {
152  sigma = segmentPathOvc_[iSeg] * segmentSigmaMom_[iSeg] /
153  (segmentMom2_[iSeg] * sqrt(segmentMom2_[iSeg] + 1 / mass_inv2) * mass_inv2);
154  sigmaTofs_.push_back(sigma);
155 
156  sigmatof += sigma * sigma;
157  }
158 
159  // compute sigma on sum of tofs assuming full correlation between segments
160  for (uint32_t iSeg = 0; iSeg < nSegment_; iSeg++) {
161  for (uint32_t jSeg = iSeg + 1; jSeg < nSegment_; jSeg++) {
162  sigmatof += 2 * sigmaTofs_[iSeg] * sigmaTofs_[jSeg];
163  }
164  }
165 
166  return sqrt(sigmatof);
167  }
168 
169  inline uint32_t size() const { return nSegment_; }
170 
171  inline uint32_t removeFirstSegment() {
172  if (nSegment_ > 0) {
173  segmentPathOvc_.erase(segmentPathOvc_.begin());
174  segmentMom2_.erase(segmentMom2_.begin());
175  nSegment_--;
176  }
177  return nSegment_;
178  }
179 
180  inline std::pair<float, float> getSegmentPathAndMom2(uint32_t iSegment) const {
181  if (iSegment >= nSegment_) {
182  throw cms::Exception("TrackExtenderWithMTD") << "Requesting non existing track segment #" << iSegment;
183  }
184  return std::make_pair(segmentPathOvc_[iSegment], segmentMom2_[iSegment]);
185  }
186 
187  uint32_t nSegment_ = 0;
188  std::vector<float> segmentPathOvc_;
189  std::vector<float> segmentMom2_;
190  std::vector<float> segmentSigmaMom_;
191 
192  std::vector<float> sigmaTofs_;
193  };
194 
195  struct TrackTofPidInfo {
196  float tmtd;
197  float tmtderror;
198  float pathlength;
199 
200  float betaerror;
201 
202  float dt;
203  float dterror;
204  float dterror2;
205  float dtchi2;
206 
207  float dt_best;
208  float dterror_best;
209  float dtchi2_best;
210 
211  float gammasq_pi;
212  float beta_pi;
213  float dt_pi;
214  float sigma_dt_pi;
215 
216  float gammasq_k;
217  float beta_k;
218  float dt_k;
219  float sigma_dt_k;
220 
221  float gammasq_p;
222  float beta_p;
223  float dt_p;
224  float sigma_dt_p;
225 
226  float prob_pi;
227  float prob_k;
228  float prob_p;
229  };
230 
231  enum class TofCalc { kCost = 1, kSegm = 2, kMixd = 3 };
232  enum class SigmaTofCalc { kCost = 1, kSegm = 2, kMixd = 3 };
233 
234  const TrackTofPidInfo computeTrackTofPidInfo(float magp2,
235  float length,
236  TrackSegments trs,
237  float t_mtd,
238  float t_mtderr,
239  float t_vtx,
240  float t_vtx_err,
241  bool addPIDError = true,
242  TofCalc choice = TofCalc::kCost,
243  SigmaTofCalc sigma_choice = SigmaTofCalc::kCost) {
244  constexpr float m_pi = 0.13957018f;
245  constexpr float m_pi_inv2 = 1.0f / m_pi / m_pi;
246  constexpr float m_k = 0.493677f;
247  constexpr float m_k_inv2 = 1.0f / m_k / m_k;
248  constexpr float m_p = 0.9382720813f;
249  constexpr float m_p_inv2 = 1.0f / m_p / m_p;
250 
251  TrackTofPidInfo tofpid;
252 
253  tofpid.tmtd = t_mtd;
254  tofpid.tmtderror = t_mtderr;
255  tofpid.pathlength = length;
256 
257  auto deltat = [&](const float mass_inv2, const float betatmp) {
258  float res(1.f);
259  switch (choice) {
260  case TofCalc::kCost:
261  res = tofpid.pathlength / betatmp * c_inv;
262  break;
263  case TofCalc::kSegm:
264  res = trs.computeTof(mass_inv2);
265  break;
266  case TofCalc::kMixd:
267  res = trs.computeTof(mass_inv2) + tofpid.pathlength / betatmp * c_inv;
268  break;
269  }
270  return res;
271  };
272 
273  auto sigmadeltat = [&](const float mass_inv2) {
274  float res(1.f);
275  switch (sigma_choice) {
276  case SigmaTofCalc::kCost:
277  // sigma(t) = sigma(p) * |dt/dp| = sigma(p) * DeltaL/c * m^2 / (p^2 * E)
278  res = tofpid.pathlength * c_inv * trs.segmentSigmaMom_[trs.nSegment_ - 1] /
279  (magp2 * sqrt(magp2 + 1 / mass_inv2) * mass_inv2);
280  break;
281  case SigmaTofCalc::kSegm:
282  res = trs.computeSigmaTof(mass_inv2);
283  break;
284  case SigmaTofCalc::kMixd:
285  float res1 = tofpid.pathlength * c_inv * trs.segmentSigmaMom_[trs.nSegment_ - 1] /
286  (magp2 * sqrt(magp2 + 1 / mass_inv2) * mass_inv2);
287  float res2 = trs.computeSigmaTof(mass_inv2);
288  res = sqrt(res1 * res1 + res2 * res2 + 2 * res1 * res2);
289  }
290 
291  return res;
292  };
293 
294  tofpid.gammasq_pi = 1.f + magp2 * m_pi_inv2;
295  tofpid.beta_pi = std::sqrt(1.f - 1.f / tofpid.gammasq_pi);
296  tofpid.dt_pi = deltat(m_pi_inv2, tofpid.beta_pi);
297  tofpid.sigma_dt_pi = sigmadeltat(m_pi_inv2);
298 
299  tofpid.gammasq_k = 1.f + magp2 * m_k_inv2;
300  tofpid.beta_k = std::sqrt(1.f - 1.f / tofpid.gammasq_k);
301  tofpid.dt_k = deltat(m_k_inv2, tofpid.beta_k);
302  tofpid.sigma_dt_k = sigmadeltat(m_k_inv2);
303 
304  tofpid.gammasq_p = 1.f + magp2 * m_p_inv2;
305  tofpid.beta_p = std::sqrt(1.f - 1.f / tofpid.gammasq_p);
306  tofpid.dt_p = deltat(m_p_inv2, tofpid.beta_p);
307  tofpid.sigma_dt_p = sigmadeltat(m_p_inv2);
308 
309  tofpid.dt = tofpid.tmtd - tofpid.dt_pi - t_vtx; //assume by default the pi hypothesis
310  tofpid.dterror2 = tofpid.tmtderror * tofpid.tmtderror + t_vtx_err * t_vtx_err;
311  tofpid.betaerror = 0.f;
312  if (addPIDError) {
313  tofpid.dterror2 = tofpid.dterror2 + (tofpid.dt_p - tofpid.dt_pi) * (tofpid.dt_p - tofpid.dt_pi);
314  tofpid.betaerror = tofpid.beta_p - tofpid.beta_pi;
315  } else {
316  // only add sigma(TOF) if not considering mass hp. uncertainty
317  tofpid.dterror2 = tofpid.dterror2 + tofpid.sigma_dt_pi * tofpid.sigma_dt_pi;
318  }
319  tofpid.dterror = sqrt(tofpid.dterror2);
320 
321  tofpid.dtchi2 = (tofpid.dt * tofpid.dt) / tofpid.dterror2;
322 
323  tofpid.dt_best = tofpid.dt;
324  tofpid.dterror_best = tofpid.dterror;
325  tofpid.dtchi2_best = tofpid.dtchi2;
326 
327  tofpid.prob_pi = -1.f;
328  tofpid.prob_k = -1.f;
329  tofpid.prob_p = -1.f;
330 
331  if (!addPIDError) {
332  //*TODO* deal with heavier nucleons and/or BSM case here?
333  const float dterror2_wo_sigmatof = tofpid.dterror2 - tofpid.sigma_dt_pi * tofpid.sigma_dt_pi;
334  float chi2_pi = tofpid.dtchi2;
335  float chi2_k = (tofpid.tmtd - tofpid.dt_k - t_vtx) * (tofpid.tmtd - tofpid.dt_k - t_vtx) /
336  (dterror2_wo_sigmatof + tofpid.sigma_dt_k * tofpid.sigma_dt_k);
337  float chi2_p = (tofpid.tmtd - tofpid.dt_p - t_vtx) * (tofpid.tmtd - tofpid.dt_p - t_vtx) /
338  (dterror2_wo_sigmatof + tofpid.sigma_dt_p * tofpid.sigma_dt_p);
339 
340  float rawprob_pi = exp(-0.5f * chi2_pi);
341  float rawprob_k = exp(-0.5f * chi2_k);
342  float rawprob_p = exp(-0.5f * chi2_p);
343  float normprob = 1.f / (rawprob_pi + rawprob_k + rawprob_p);
344 
345  tofpid.prob_pi = rawprob_pi * normprob;
346  tofpid.prob_k = rawprob_k * normprob;
347  tofpid.prob_p = rawprob_p * normprob;
348 
349  float prob_heavy = 1.f - tofpid.prob_pi;
350  constexpr float heavy_threshold = 0.75f;
351 
352  if (prob_heavy > heavy_threshold) {
353  if (chi2_k < chi2_p) {
354  tofpid.dt_best = (tofpid.tmtd - tofpid.dt_k - t_vtx);
355  tofpid.dtchi2_best = chi2_k;
356  } else {
357  tofpid.dt_best = (tofpid.tmtd - tofpid.dt_p - t_vtx);
358  tofpid.dtchi2_best = chi2_p;
359  }
360  }
361  }
362  return tofpid;
363  }
364 
365  bool getTrajectoryStateClosestToBeamLine(const Trajectory& traj,
366  const reco::BeamSpot& bs,
367  const Propagator* thePropagator,
369  // get the state closest to the beamline
370  TrajectoryStateOnSurface stateForProjectionToBeamLineOnSurface =
371  traj.closestMeasurement(GlobalPoint(bs.x0(), bs.y0(), bs.z0())).updatedState();
372 
373  if (!stateForProjectionToBeamLineOnSurface.isValid()) {
374  edm::LogError("CannotPropagateToBeamLine") << "the state on the closest measurement isnot valid. skipping track.";
375  return false;
376  }
377 
378  const FreeTrajectoryState& stateForProjectionToBeamLine = *stateForProjectionToBeamLineOnSurface.freeState();
379 
380  TSCBLBuilderWithPropagator tscblBuilder(*thePropagator);
381  tscbl = tscblBuilder(stateForProjectionToBeamLine, bs);
382 
383  return tscbl.isValid();
384  }
385 
386  bool trackPathLength(const Trajectory& traj,
388  const Propagator* thePropagator,
389  float& pathlength,
390  TrackSegments& trs) {
391  pathlength = 0.f;
392 
393  bool validpropagation = true;
394  float oldp = traj.measurements().begin()->updatedState().globalMomentum().mag();
395  float pathlength1 = 0.f;
396  float pathlength2 = 0.f;
397 
398  //add pathlength layer by layer
399  for (auto it = traj.measurements().begin(); it != traj.measurements().end() - 1; ++it) {
400  const auto& propresult = thePropagator->propagateWithPath(it->updatedState(), (it + 1)->updatedState().surface());
401  float layerpathlength = std::abs(propresult.second);
402  if (layerpathlength == 0.f) {
403  validpropagation = false;
404  }
405  pathlength1 += layerpathlength;
406 
407  // sigma(p) from curvilinear error (on q/p)
408  float sigma_p = sqrt((it + 1)->updatedState().curvilinearError().matrix()(0, 0)) *
409  (it + 1)->updatedState().globalMomentum().mag2();
410 
411  trs.addSegment(layerpathlength, (it + 1)->updatedState().globalMomentum().mag2(), sigma_p);
412 
413  LogTrace("TrackExtenderWithMTD") << "TSOS " << std::fixed << std::setw(4) << trs.size() << " R_i " << std::fixed
414  << std::setw(14) << it->updatedState().globalPosition().perp() << " z_i "
415  << std::fixed << std::setw(14) << it->updatedState().globalPosition().z()
416  << " R_e " << std::fixed << std::setw(14)
417  << (it + 1)->updatedState().globalPosition().perp() << " z_e " << std::fixed
418  << std::setw(14) << (it + 1)->updatedState().globalPosition().z() << " p "
419  << std::fixed << std::setw(14) << (it + 1)->updatedState().globalMomentum().mag()
420  << " dp " << std::fixed << std::setw(14)
421  << (it + 1)->updatedState().globalMomentum().mag() - oldp;
422  oldp = (it + 1)->updatedState().globalMomentum().mag();
423  }
424 
425  //add distance from bs to first measurement
426  auto const& tscblPCA = tscbl.trackStateAtPCA();
427  auto const& aSurface = traj.direction() == alongMomentum ? traj.firstMeasurement().updatedState().surface()
428  : traj.lastMeasurement().updatedState().surface();
429  pathlength2 = thePropagator->propagateWithPath(tscblPCA, aSurface).second;
430  if (pathlength2 == 0.f) {
431  validpropagation = false;
432  }
433  pathlength = pathlength1 + pathlength2;
434 
435  float sigma_p = sqrt(tscblPCA.curvilinearError().matrix()(0, 0)) * tscblPCA.momentum().mag2();
436 
437  trs.addSegment(pathlength2, tscblPCA.momentum().mag2(), sigma_p);
438 
439  LogTrace("TrackExtenderWithMTD") << "TSOS " << std::fixed << std::setw(4) << trs.size() << " R_e " << std::fixed
440  << std::setw(14) << tscblPCA.position().perp() << " z_e " << std::fixed
441  << std::setw(14) << tscblPCA.position().z() << " p " << std::fixed << std::setw(14)
442  << tscblPCA.momentum().mag() << " dp " << std::fixed << std::setw(14)
443  << tscblPCA.momentum().mag() - oldp << " sigma_p = " << std::fixed << std::setw(14)
444  << sigma_p << " sigma_p/p = " << std::fixed << std::setw(14)
445  << sigma_p / tscblPCA.momentum().mag() * 100 << " %";
446 
447  return validpropagation;
448  }
449 
450  bool trackPathLength(const Trajectory& traj,
451  const reco::BeamSpot& bs,
452  const Propagator* thePropagator,
453  float& pathlength,
454  TrackSegments& trs) {
455  pathlength = 0.f;
456 
458  bool tscbl_status = getTrajectoryStateClosestToBeamLine(traj, bs, thePropagator, tscbl);
459 
460  if (!tscbl_status)
461  return false;
462 
463  return trackPathLength(traj, tscbl, thePropagator, pathlength, trs);
464  }
465 
466 } // namespace
467 
468 template <class TrackCollection>
470 public:
473 
475 
476  template <class H, class T>
477  void fillValueMap(edm::Event& iEvent, const H& handle, const std::vector<T>& vec, const edm::EDPutToken& token) const;
478 
479  void produce(edm::Event& ev, const edm::EventSetup& es) final;
480 
481  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
482 
484  const Trajectory& traj,
485  const float,
486  const float,
487  const TrackSegments&,
489  const MTDDetLayerGeometry*,
490  const MagneticField* field,
491  const Propagator* prop,
492  const reco::BeamSpot& bs,
493  const float vtxTime,
494  const bool matchVertex,
495  MTDHitMatchingInfo& bestHit) const;
496 
498  const Trajectory& traj,
499  const float,
500  const float,
501  const TrackSegments&,
503  const MTDDetLayerGeometry*,
504  const MagneticField* field,
505  const Propagator* prop,
506  const reco::BeamSpot& bs,
507  const float vtxTime,
508  const bool matchVertex,
509  MTDHitMatchingInfo& bestHit) const;
510 
511  void fillMatchingHits(const DetLayer*,
513  const Trajectory&,
514  const float,
515  const float,
516  const TrackSegments&,
518  const Propagator*,
519  const reco::BeamSpot&,
520  const float&,
521  const bool,
523  MTDHitMatchingInfo&) const;
524 
527  if (!recHits.empty()) {
528  GlobalPoint first = gtg_->idToDet(recHits.front()->geographicalId())->position();
529  GlobalPoint last = gtg_->idToDet(recHits.back()->geographicalId())->position();
530 
531  // maybe perp2?
532  auto rFirst = first.mag2();
533  auto rLast = last.mag2();
534  if (rFirst < rLast)
536  if (rFirst > rLast)
538  }
539  LogDebug("TrackExtenderWithMTD") << "Impossible to determine the rechits order" << endl;
541  }
542 
543  reco::Track buildTrack(const reco::TrackRef&,
544  const Trajectory&,
545  const Trajectory&,
546  const reco::BeamSpot&,
547  const MagneticField* field,
548  const Propagator* prop,
549  bool hasMTD,
550  float& pathLength,
551  float& tmtdOut,
552  float& sigmatmtdOut,
553  GlobalPoint& tmtdPosOut,
554  float& tofpi,
555  float& tofk,
556  float& tofp,
557  float& sigmatofpi,
558  float& sigmatofk,
559  float& sigmatofp) const;
560  reco::TrackExtra buildTrackExtra(const Trajectory& trajectory) const;
561 
562  string dumpLayer(const DetLayer* layer) const;
563 
564 private:
587 
595 
596  const bool updateTraj_, updateExtra_, updatePattern_;
597  const std::string mtdRecHitBuilder_, propagator_, transientTrackBuilder_;
598  std::unique_ptr<MeasurementEstimator> theEstimator;
599  std::unique_ptr<TrackTransformer> theTransformer;
606 
611 
612  const float estMaxChi2_;
613  const float estMaxNSigma_;
614  const float btlChi2Cut_;
615  const float btlTimeChi2Cut_;
616  const float etlChi2Cut_;
617  const float etlTimeChi2Cut_;
618 
619  const bool useVertex_;
620  const bool useSimVertex_;
621  const float dzCut_;
622  const float bsTimeSpread_;
623 };
624 
625 template <class TrackCollection>
627  : tracksToken_(consumes<InputCollection>(iConfig.getParameter<edm::InputTag>("tracksSrc"))),
628  trajTrackAToken_(consumes<TrajTrackAssociationCollection>(iConfig.getParameter<edm::InputTag>("trjtrkAssSrc"))),
629  hitsToken_(consumes<MTDTrackingDetSetVector>(iConfig.getParameter<edm::InputTag>("hitsSrc"))),
630  bsToken_(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpotSrc"))),
631  updateTraj_(iConfig.getParameter<bool>("updateTrackTrajectory")),
632  updateExtra_(iConfig.getParameter<bool>("updateTrackExtra")),
633  updatePattern_(iConfig.getParameter<bool>("updateTrackHitPattern")),
634  mtdRecHitBuilder_(iConfig.getParameter<std::string>("MTDRecHitBuilder")),
635  propagator_(iConfig.getParameter<std::string>("Propagator")),
636  transientTrackBuilder_(iConfig.getParameter<std::string>("TransientTrackBuilder")),
637  estMaxChi2_(iConfig.getParameter<double>("estimatorMaxChi2")),
638  estMaxNSigma_(iConfig.getParameter<double>("estimatorMaxNSigma")),
639  btlChi2Cut_(iConfig.getParameter<double>("btlChi2Cut")),
640  btlTimeChi2Cut_(iConfig.getParameter<double>("btlTimeChi2Cut")),
641  etlChi2Cut_(iConfig.getParameter<double>("etlChi2Cut")),
642  etlTimeChi2Cut_(iConfig.getParameter<double>("etlTimeChi2Cut")),
643  useVertex_(iConfig.getParameter<bool>("useVertex")),
644  useSimVertex_(iConfig.getParameter<bool>("useSimVertex")),
645  dzCut_(iConfig.getParameter<double>("dZCut")),
646  bsTimeSpread_(iConfig.getParameter<double>("bsTimeSpread")) {
647  if (useVertex_) {
648  if (useSimVertex_) {
649  genVtxPositionToken_ = consumes<GlobalPoint>(iConfig.getParameter<edm::InputTag>("genVtxPositionSrc"));
650  genVtxTimeToken_ = consumes<float>(iConfig.getParameter<edm::InputTag>("genVtxTimeSrc"));
651  } else
652  vtxToken_ = consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vtxSrc"));
653  }
654 
655  theEstimator = std::make_unique<Chi2MeasurementEstimator>(estMaxChi2_, estMaxNSigma_);
656  theTransformer = std::make_unique<TrackTransformer>(iConfig.getParameterSet("TrackTransformer"), consumesCollector());
657 
658  btlMatchChi2Token_ = produces<edm::ValueMap<float>>("btlMatchChi2");
659  etlMatchChi2Token_ = produces<edm::ValueMap<float>>("etlMatchChi2");
660  btlMatchTimeChi2Token_ = produces<edm::ValueMap<float>>("btlMatchTimeChi2");
661  etlMatchTimeChi2Token_ = produces<edm::ValueMap<float>>("etlMatchTimeChi2");
662  npixBarrelToken_ = produces<edm::ValueMap<int>>("npixBarrel");
663  npixEndcapToken_ = produces<edm::ValueMap<int>>("npixEndcap");
664  outermostHitPositionToken_ = produces<edm::ValueMap<float>>("generalTrackOutermostHitPosition");
665  pOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackp");
666  betaOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackBeta");
667  t0OrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackt0");
668  sigmat0OrigTrkToken_ = produces<edm::ValueMap<float>>("generalTracksigmat0");
669  pathLengthOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackPathLength");
670  tmtdOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTracktmtd");
671  sigmatmtdOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTracksigmatmtd");
672  tmtdPosOrigTrkToken_ = produces<edm::ValueMap<GlobalPoint>>("generalTrackmtdpos");
673  tofpiOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackTofPi");
674  tofkOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackTofK");
675  tofpOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackTofP");
676  sigmatofpiOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackSigmaTofPi");
677  sigmatofkOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackSigmaTofK");
678  sigmatofpOrigTrkToken_ = produces<edm::ValueMap<float>>("generalTrackSigmaTofP");
679  assocOrigTrkToken_ = produces<edm::ValueMap<int>>("generalTrackassoc");
680 
681  builderToken_ = esConsumes<TransientTrackBuilder, TransientTrackRecord>(edm::ESInputTag("", transientTrackBuilder_));
683  esConsumes<TransientTrackingRecHitBuilder, TransientRecHitRecord>(edm::ESInputTag("", mtdRecHitBuilder_));
684  gtgToken_ = esConsumes<GlobalTrackingGeometry, GlobalTrackingGeometryRecord>();
685  dlgeoToken_ = esConsumes<MTDDetLayerGeometry, MTDRecoGeometryRecord>();
686  magfldToken_ = esConsumes<MagneticField, IdealMagneticFieldRecord>();
687  propToken_ = esConsumes<Propagator, TrackingComponentsRecord>(edm::ESInputTag("", propagator_));
688  ttopoToken_ = esConsumes<TrackerTopology, TrackerTopologyRcd>();
689 
690  produces<edm::OwnVector<TrackingRecHit>>();
691  produces<reco::TrackExtraCollection>();
692  produces<TrackCollection>();
693 }
694 
695 template <class TrackCollection>
698  desc.add<edm::InputTag>("tracksSrc", edm::InputTag("generalTracks"));
699  desc.add<edm::InputTag>("trjtrkAssSrc", edm::InputTag("generalTracks"));
700  desc.add<edm::InputTag>("hitsSrc", edm::InputTag("mtdTrackingRecHits"));
701  desc.add<edm::InputTag>("beamSpotSrc", edm::InputTag("offlineBeamSpot"));
702  desc.add<edm::InputTag>("genVtxPositionSrc", edm::InputTag("genParticles:xyz0"));
703  desc.add<edm::InputTag>("genVtxTimeSrc", edm::InputTag("genParticles:t0"));
704  desc.add<edm::InputTag>("vtxSrc", edm::InputTag("offlinePrimaryVertices4D"));
705  desc.add<bool>("updateTrackTrajectory", true);
706  desc.add<bool>("updateTrackExtra", true);
707  desc.add<bool>("updateTrackHitPattern", true);
708  desc.add<std::string>("TransientTrackBuilder", "TransientTrackBuilder");
709  desc.add<std::string>("MTDRecHitBuilder", "MTDRecHitBuilder");
710  desc.add<std::string>("Propagator", "PropagatorWithMaterialForMTD");
712  false,
713  "KFFitterForRefitInsideOut",
714  "KFSmootherForRefitInsideOut",
715  "PropagatorWithMaterialForMTD",
716  "alongMomentum",
717  true,
718  "WithTrackAngle",
719  "MuonRecHitBuilder",
720  "MTDRecHitBuilder");
721  desc.add<edm::ParameterSetDescription>("TrackTransformer", transDesc);
722  desc.add<double>("estimatorMaxChi2", 500.);
723  desc.add<double>("estimatorMaxNSigma", 10.);
724  desc.add<double>("btlChi2Cut", 50.);
725  desc.add<double>("btlTimeChi2Cut", 10.);
726  desc.add<double>("etlChi2Cut", 50.);
727  desc.add<double>("etlTimeChi2Cut", 10.);
728  desc.add<bool>("useVertex", false);
729  desc.add<bool>("useSimVertex", false);
730  desc.add<double>("dZCut", 0.1);
731  desc.add<double>("bsTimeSpread", 0.2);
732  descriptions.add("trackExtenderWithMTDBase", desc);
733 }
734 
735 template <class TrackCollection>
736 template <class H, class T>
738  const H& handle,
739  const std::vector<T>& vec,
740  const edm::EDPutToken& token) const {
741  auto out = std::make_unique<edm::ValueMap<T>>();
743  filler.insert(handle, vec.begin(), vec.end());
744  filler.fill();
745  iEvent.put(token, std::move(out));
746 }
747 
748 template <class TrackCollection>
750  //this produces pieces of the track extra
751  Traj2TrackHits t2t;
752 
753  theTransformer->setServices(es);
754  TrackingRecHitRefProd hitsRefProd = ev.getRefBeforePut<TrackingRecHitCollection>();
755  reco::TrackExtraRefProd extrasRefProd = ev.getRefBeforePut<reco::TrackExtraCollection>();
756 
757  gtg_ = es.getHandle(gtgToken_);
758 
759  auto geo = es.getTransientHandle(dlgeoToken_);
760 
761  auto magfield = es.getTransientHandle(magfldToken_);
762 
763  builder_ = es.getHandle(builderToken_);
764  hitbuilder_ = es.getHandle(hitbuilderToken_);
765 
766  auto propH = es.getTransientHandle(propToken_);
767  const Propagator* prop = propH.product();
768 
769  auto httopo = es.getTransientHandle(ttopoToken_);
770  const TrackerTopology& ttopo = *httopo;
771 
772  auto output = std::make_unique<TrackCollection>();
773  auto extras = std::make_unique<reco::TrackExtraCollection>();
774  auto outhits = std::make_unique<edm::OwnVector<TrackingRecHit>>();
775 
776  std::vector<float> btlMatchChi2;
777  std::vector<float> etlMatchChi2;
778  std::vector<float> btlMatchTimeChi2;
779  std::vector<float> etlMatchTimeChi2;
780  std::vector<int> npixBarrel;
781  std::vector<int> npixEndcap;
782  std::vector<float> outermostHitPosition;
783  std::vector<float> pOrigTrkRaw;
784  std::vector<float> betaOrigTrkRaw;
785  std::vector<float> t0OrigTrkRaw;
786  std::vector<float> sigmat0OrigTrkRaw;
787  std::vector<float> pathLengthsOrigTrkRaw;
788  std::vector<float> tmtdOrigTrkRaw;
789  std::vector<float> sigmatmtdOrigTrkRaw;
790  std::vector<GlobalPoint> tmtdPosOrigTrkRaw;
791  std::vector<float> tofpiOrigTrkRaw;
792  std::vector<float> tofkOrigTrkRaw;
793  std::vector<float> tofpOrigTrkRaw;
794  std::vector<float> sigmatofpiOrigTrkRaw;
795  std::vector<float> sigmatofkOrigTrkRaw;
796  std::vector<float> sigmatofpOrigTrkRaw;
797  std::vector<int> assocOrigTrkRaw;
798 
799  auto const tracksH = ev.getHandle(tracksToken_);
800 
801  const auto& trjtrks = ev.get(trajTrackAToken_);
802 
803  //MTD hits DetSet
804  const auto& hits = ev.get(hitsToken_);
805 
806  //beam spot
807  const auto& bs = ev.get(bsToken_);
808 
809  const Vertex* pv = nullptr;
810  if (useVertex_ && !useSimVertex_) {
811  auto const& vtxs = ev.get(vtxToken_);
812  if (!vtxs.empty())
813  pv = &vtxs[0];
814  }
815 
816  std::unique_ptr<math::XYZTLorentzVectorF> genPV(nullptr);
817  if (useVertex_ && useSimVertex_) {
818  const auto& genVtxPosition = ev.get(genVtxPositionToken_);
819  const auto& genVtxTime = ev.get(genVtxTimeToken_);
820  genPV = std::make_unique<math::XYZTLorentzVectorF>(
821  genVtxPosition.x(), genVtxPosition.y(), genVtxPosition.z(), genVtxTime);
822  }
823 
824  float vtxTime = 0.f;
825  if (useVertex_) {
826  if (useSimVertex_ && genPV) {
827  vtxTime = genPV->t();
828  } else if (pv)
829  vtxTime = pv->t(); //already in ns
830  }
831 
832  std::vector<unsigned> track_indices;
833  unsigned itrack = 0;
834 
835  for (const auto& trjtrk : trjtrks) {
836  const Trajectory& trajs = *trjtrk.key;
837  const reco::TrackRef& track = trjtrk.val;
838 
839  LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: extrapolating track " << itrack
840  << " p/pT = " << track->p() << " " << track->pt() << " eta = " << track->eta();
841  LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: sigma_p = "
842  << sqrt(track->covariance()(0, 0)) * track->p2()
843  << " sigma_p/p = " << sqrt(track->covariance()(0, 0)) * track->p() * 100 << " %";
844 
845  float trackVtxTime = 0.f;
846  if (useVertex_) {
847  float dz;
848  if (useSimVertex_)
849  dz = std::abs(track->dz(math::XYZPoint(*genPV)));
850  else
851  dz = std::abs(track->dz(pv->position()));
852 
853  if (dz < dzCut_)
854  trackVtxTime = vtxTime;
855  }
856 
857  reco::TransientTrack ttrack(track, magfield.product(), gtg_);
858  auto thits = theTransformer->getTransientRecHits(ttrack);
860  MTDHitMatchingInfo mBTL, mETL;
861 
862  if (trajs.isValid()) {
863  // get the outermost trajectory point on the track
864  TrajectoryStateOnSurface tsos = builder_->build(track).outermostMeasurementState();
866  bool tscbl_status = getTrajectoryStateClosestToBeamLine(trajs, bs, prop, tscbl);
867 
868  if (tscbl_status) {
869  float pmag2 = tscbl.trackStateAtPCA().momentum().mag2();
870  float pathlength0;
871  TrackSegments trs0;
872  trackPathLength(trajs, tscbl, prop, pathlength0, trs0);
873 
874  const auto& btlhits = tryBTLLayers(tsos,
875  trajs,
876  pmag2,
877  pathlength0,
878  trs0,
879  hits,
880  geo.product(),
881  magfield.product(),
882  prop,
883  bs,
884  trackVtxTime,
885  trackVtxTime != 0.f,
886  mBTL);
887  mtdthits.insert(mtdthits.end(), btlhits.begin(), btlhits.end());
888 
889  // in the future this should include an intermediate refit before propagating to the ETL
890  // for now it is ok
891  const auto& etlhits = tryETLLayers(tsos,
892  trajs,
893  pmag2,
894  pathlength0,
895  trs0,
896  hits,
897  geo.product(),
898  magfield.product(),
899  prop,
900  bs,
901  trackVtxTime,
902  trackVtxTime != 0.f,
903  mETL);
904  mtdthits.insert(mtdthits.end(), etlhits.begin(), etlhits.end());
905  }
906 #ifdef EDM_ML_DEBUG
907  else {
908  LogTrace("TrackExtenderWithMTD") << "Failing getTrajectoryStateClosestToBeamLine, no search for hits in MTD!";
909  }
910 #endif
911  }
912 
913  auto ordering = checkRecHitsOrdering(thits);
915  thits.insert(thits.end(), mtdthits.begin(), mtdthits.end());
916  } else {
917  std::reverse(mtdthits.begin(), mtdthits.end());
918  mtdthits.insert(mtdthits.end(), thits.begin(), thits.end());
919  thits.swap(mtdthits);
920  }
921 
922  const auto& trajwithmtd =
923  mtdthits.empty() ? std::vector<Trajectory>(1, trajs) : theTransformer->transform(ttrack, thits);
924  float pMap = 0.f, betaMap = 0.f, t0Map = 0.f, sigmat0Map = -1.f, pathLengthMap = -1.f, tmtdMap = 0.f,
925  sigmatmtdMap = -1.f, tofpiMap = 0.f, tofkMap = 0.f, tofpMap = 0.f, sigmatofpiMap = -1.f, sigmatofkMap = -1.f,
926  sigmatofpMap = -1.f;
927  GlobalPoint tmtdPosMap{0., 0., 0.};
928  int iMap = -1;
929 
930  for (const auto& trj : trajwithmtd) {
931  const auto& thetrj = (updateTraj_ ? trj : trajs);
932  float pathLength = 0.f, tmtd = 0.f, sigmatmtd = -1.f, tofpi = 0.f, tofk = 0.f, tofp = 0.f, sigmatofpi = -1.f,
933  sigmatofk = -1.f, sigmatofp = -1.f;
934  GlobalPoint tmtdPos{0., 0., 0.};
935  LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: refit track " << itrack << " p/pT = " << track->p()
936  << " " << track->pt() << " eta = " << track->eta();
937  reco::Track result = buildTrack(track,
938  thetrj,
939  trj,
940  bs,
941  magfield.product(),
942  prop,
943  !trajwithmtd.empty() && !mtdthits.empty(),
944  pathLength,
945  tmtd,
946  sigmatmtd,
947  tmtdPos,
948  tofpi,
949  tofk,
950  tofp,
951  sigmatofpi,
952  sigmatofk,
953  sigmatofp);
954  if (result.ndof() >= 0) {
956  reco::TrackExtra::TrajParams trajParams;
958  size_t hitsstart = outhits->size();
959  if (updatePattern_) {
960  t2t(trj, *outhits, trajParams, chi2s); // this fills the output hit collection
961  } else {
962  t2t(thetrj, *outhits, trajParams, chi2s);
963  }
964  size_t hitsend = outhits->size();
965  extras->push_back(buildTrackExtra(trj)); // always push back the fully built extra, update by setting in track
966  extras->back().setHits(hitsRefProd, hitsstart, hitsend - hitsstart);
967  extras->back().setTrajParams(trajParams, chi2s);
968  //create the track
969  output->push_back(result);
970  btlMatchChi2.push_back(mBTL.hit ? mBTL.estChi2 : -1.f);
971  etlMatchChi2.push_back(mETL.hit ? mETL.estChi2 : -1.f);
972  btlMatchTimeChi2.push_back(mBTL.hit ? mBTL.timeChi2 : -1.f);
973  etlMatchTimeChi2.push_back(mETL.hit ? mETL.timeChi2 : -1.f);
974  pathLengthMap = pathLength;
975  tmtdMap = tmtd;
976  sigmatmtdMap = sigmatmtd;
977  tmtdPosMap = tmtdPos;
978  auto& backtrack = output->back();
979  iMap = output->size() - 1;
980  pMap = backtrack.p();
981  betaMap = backtrack.beta();
982  t0Map = backtrack.t0();
983  sigmat0Map = std::copysign(std::sqrt(std::abs(backtrack.covt0t0())), backtrack.covt0t0());
984  tofpiMap = tofpi;
985  tofkMap = tofk;
986  tofpMap = tofp;
987  sigmatofpiMap = sigmatofpi;
988  sigmatofkMap = sigmatofk;
989  sigmatofpMap = sigmatofp;
990  reco::TrackExtraRef extraRef(extrasRefProd, extras->size() - 1);
991  backtrack.setExtra((updateExtra_ ? extraRef : track->extra()));
992  for (unsigned ihit = hitsstart; ihit < hitsend; ++ihit) {
993  backtrack.appendHitPattern((*outhits)[ihit], ttopo);
994  }
995 #ifdef EDM_ML_DEBUG
996  LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: hit pattern of refitted track";
997  for (int i = 0; i < backtrack.hitPattern().numberOfAllHits(reco::HitPattern::TRACK_HITS); i++) {
998  backtrack.hitPattern().printHitPattern(reco::HitPattern::TRACK_HITS, i, std::cout);
999  }
1000  LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: missing hit pattern of refitted track";
1001  for (int i = 0; i < backtrack.hitPattern().numberOfAllHits(reco::HitPattern::MISSING_INNER_HITS); i++) {
1002  backtrack.hitPattern().printHitPattern(reco::HitPattern::MISSING_INNER_HITS, i, std::cout);
1003  }
1004 #endif
1005  npixBarrel.push_back(backtrack.hitPattern().numberOfValidPixelBarrelHits());
1006  npixEndcap.push_back(backtrack.hitPattern().numberOfValidPixelEndcapHits());
1007  outermostHitPosition.push_back(
1008  mBTL.hit ? (float)(*track).outerRadius()
1009  : (float)(*track).outerZ()); // save R of the outermost hit for BTL, z for ETL.
1010 
1011  LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: tmtd " << tmtdMap << " +/- " << sigmatmtdMap
1012  << " t0 " << t0Map << " +/- " << sigmat0Map << " tof pi/K/p " << tofpiMap
1013  << "+/-" << fmt::format("{:0.2g}", sigmatofpiMap) << " ("
1014  << fmt::format("{:0.2g}", sigmatofpiMap / tofpiMap * 100) << "%) " << tofkMap
1015  << "+/-" << fmt::format("{:0.2g}", sigmatofkMap) << " ("
1016  << fmt::format("{:0.2g}", sigmatofkMap / tofkMap * 100) << "%) " << tofpMap
1017  << "+/-" << fmt::format("{:0.2g}", sigmatofpMap) << " ("
1018  << fmt::format("{:0.2g}", sigmatofpMap / tofpMap * 100) << "%) ";
1019  } else {
1020  LogTrace("TrackExtenderWithMTD") << "Error in the MTD track refitting. This should not happen";
1021  }
1022  }
1023 
1024  pOrigTrkRaw.push_back(pMap);
1025  betaOrigTrkRaw.push_back(betaMap);
1026  t0OrigTrkRaw.push_back(t0Map);
1027  sigmat0OrigTrkRaw.push_back(sigmat0Map);
1028  pathLengthsOrigTrkRaw.push_back(pathLengthMap);
1029  tmtdOrigTrkRaw.push_back(tmtdMap);
1030  sigmatmtdOrigTrkRaw.push_back(sigmatmtdMap);
1031  tmtdPosOrigTrkRaw.push_back(tmtdPosMap);
1032  tofpiOrigTrkRaw.push_back(tofpiMap);
1033  tofkOrigTrkRaw.push_back(tofkMap);
1034  tofpOrigTrkRaw.push_back(tofpMap);
1035  sigmatofpiOrigTrkRaw.push_back(sigmatofpiMap);
1036  sigmatofkOrigTrkRaw.push_back(sigmatofkMap);
1037  sigmatofpOrigTrkRaw.push_back(sigmatofpMap);
1038  assocOrigTrkRaw.push_back(iMap);
1039 
1040  if (iMap == -1) {
1041  btlMatchChi2.push_back(-1.f);
1042  etlMatchChi2.push_back(-1.f);
1043  btlMatchTimeChi2.push_back(-1.f);
1044  etlMatchTimeChi2.push_back(-1.f);
1045  npixBarrel.push_back(-1.f);
1046  npixEndcap.push_back(-1.f);
1047  outermostHitPosition.push_back(0.);
1048  }
1049 
1050  ++itrack;
1051  }
1052 
1053  auto outTrksHandle = ev.put(std::move(output));
1054  ev.put(std::move(extras));
1055  ev.put(std::move(outhits));
1056 
1057  fillValueMap(ev, tracksH, btlMatchChi2, btlMatchChi2Token_);
1058  fillValueMap(ev, tracksH, etlMatchChi2, etlMatchChi2Token_);
1059  fillValueMap(ev, tracksH, btlMatchTimeChi2, btlMatchTimeChi2Token_);
1060  fillValueMap(ev, tracksH, etlMatchTimeChi2, etlMatchTimeChi2Token_);
1061  fillValueMap(ev, tracksH, npixBarrel, npixBarrelToken_);
1062  fillValueMap(ev, tracksH, npixEndcap, npixEndcapToken_);
1063  fillValueMap(ev, tracksH, outermostHitPosition, outermostHitPositionToken_);
1064  fillValueMap(ev, tracksH, pOrigTrkRaw, pOrigTrkToken_);
1065  fillValueMap(ev, tracksH, betaOrigTrkRaw, betaOrigTrkToken_);
1066  fillValueMap(ev, tracksH, t0OrigTrkRaw, t0OrigTrkToken_);
1067  fillValueMap(ev, tracksH, sigmat0OrigTrkRaw, sigmat0OrigTrkToken_);
1068  fillValueMap(ev, tracksH, pathLengthsOrigTrkRaw, pathLengthOrigTrkToken_);
1069  fillValueMap(ev, tracksH, tmtdOrigTrkRaw, tmtdOrigTrkToken_);
1070  fillValueMap(ev, tracksH, sigmatmtdOrigTrkRaw, sigmatmtdOrigTrkToken_);
1071  fillValueMap(ev, tracksH, tmtdPosOrigTrkRaw, tmtdPosOrigTrkToken_);
1072  fillValueMap(ev, tracksH, tofpiOrigTrkRaw, tofpiOrigTrkToken_);
1073  fillValueMap(ev, tracksH, tofkOrigTrkRaw, tofkOrigTrkToken_);
1074  fillValueMap(ev, tracksH, tofpOrigTrkRaw, tofpOrigTrkToken_);
1075  fillValueMap(ev, tracksH, sigmatofpiOrigTrkRaw, sigmatofpiOrigTrkToken_);
1076  fillValueMap(ev, tracksH, sigmatofkOrigTrkRaw, sigmatofkOrigTrkToken_);
1077  fillValueMap(ev, tracksH, sigmatofpOrigTrkRaw, sigmatofpOrigTrkToken_);
1078  fillValueMap(ev, tracksH, assocOrigTrkRaw, assocOrigTrkToken_);
1079 }
1080 
1081 namespace {
1082  bool cmp_for_detset(const unsigned one, const unsigned two) { return one < two; };
1083 
1084  void find_hits_in_dets(const MTDTrackingDetSetVector& hits,
1085  const Trajectory& traj,
1086  const DetLayer* layer,
1087  const TrajectoryStateOnSurface& tsos,
1088  const float pmag2,
1089  const float pathlength0,
1090  const TrackSegments& trs0,
1091  const float vtxTime,
1092  const reco::BeamSpot& bs,
1093  const float bsTimeSpread,
1094  const Propagator* prop,
1096  bool useVtxConstraint,
1097  std::set<MTDHitMatchingInfo>& out) {
1098  pair<bool, TrajectoryStateOnSurface> comp = layer->compatible(tsos, *prop, *estimator);
1099  if (comp.first) {
1100  const vector<DetLayer::DetWithState> compDets = layer->compatibleDets(tsos, *prop, *estimator);
1101  LogTrace("TrackExtenderWithMTD") << "Hit search: Compatible dets " << compDets.size();
1102  if (!compDets.empty()) {
1103  for (const auto& detWithState : compDets) {
1104  auto range = hits.equal_range(detWithState.first->geographicalId(), cmp_for_detset);
1105  if (range.first == range.second) {
1106  LogTrace("TrackExtenderWithMTD")
1107  << "Hit search: no hit in DetId " << detWithState.first->geographicalId().rawId();
1108  continue;
1109  }
1110 
1111  auto pl = prop->propagateWithPath(tsos, detWithState.second.surface());
1112  if (pl.second == 0.) {
1113  LogTrace("TrackExtenderWithMTD")
1114  << "Hit search: no propagation to DetId " << detWithState.first->geographicalId().rawId();
1115  continue;
1116  }
1117 
1118  const float t_vtx = useVtxConstraint ? vtxTime : 0.f;
1119 
1120  constexpr float vtx_res = 0.008f;
1121  const float t_vtx_err = useVtxConstraint ? vtx_res : bsTimeSpread;
1122 
1123  float lastpmag2 = trs0.getSegmentPathAndMom2(0).second;
1124 
1125  for (auto detitr = range.first; detitr != range.second; ++detitr) {
1126  for (const auto& hit : *detitr) {
1127  auto est = estimator->estimate(detWithState.second, hit);
1128  if (!est.first) {
1129  LogTrace("TrackExtenderWithMTD")
1130  << "Hit search: no compatible estimate in DetId " << detWithState.first->geographicalId().rawId()
1131  << " for hit at pos (" << std::fixed << std::setw(14) << hit.globalPosition().x() << ","
1132  << std::fixed << std::setw(14) << hit.globalPosition().y() << "," << std::fixed << std::setw(14)
1133  << hit.globalPosition().z() << ")";
1134  continue;
1135  }
1136 
1137  LogTrace("TrackExtenderWithMTD")
1138  << "Hit search: spatial compatibility DetId " << detWithState.first->geographicalId().rawId()
1139  << " TSOS dx/dy " << std::fixed << std::setw(14)
1140  << std::sqrt(detWithState.second.localError().positionError().xx()) << " " << std::fixed
1141  << std::setw(14) << std::sqrt(detWithState.second.localError().positionError().yy()) << " hit dx/dy "
1142  << std::fixed << std::setw(14) << std::sqrt(hit.localPositionError().xx()) << " " << std::fixed
1143  << std::setw(14) << std::sqrt(hit.localPositionError().yy()) << " chi2 " << std::fixed
1144  << std::setw(14) << est.second;
1145 
1146  TrackTofPidInfo tof = computeTrackTofPidInfo(lastpmag2,
1147  std::abs(pl.second),
1148  trs0,
1149  hit.time(),
1150  hit.timeError(),
1151  t_vtx,
1152  t_vtx_err, //put vtx error by hand for the moment
1153  false,
1154  TofCalc::kMixd,
1155  SigmaTofCalc::kMixd);
1156  MTDHitMatchingInfo mi;
1157  mi.hit = &hit;
1158  mi.estChi2 = est.second;
1159  mi.timeChi2 = tof.dtchi2_best; //use the chi2 for the best matching hypothesis
1160 
1161  out.insert(mi);
1162  }
1163  }
1164  }
1165  }
1166  }
1167  }
1168 } // namespace
1169 
1170 template <class TrackCollection>
1172  const TrajectoryStateOnSurface& tsos,
1173  const Trajectory& traj,
1174  const float pmag2,
1175  const float pathlength0,
1176  const TrackSegments& trs0,
1178  const MTDDetLayerGeometry* geo,
1179  const MagneticField* field,
1180  const Propagator* prop,
1181  const reco::BeamSpot& bs,
1182  const float vtxTime,
1183  const bool matchVertex,
1184  MTDHitMatchingInfo& bestHit) const {
1185  const vector<const DetLayer*>& layers = geo->allBTLLayers();
1186 
1188  bestHit = MTDHitMatchingInfo();
1189  for (const DetLayer* ilay : layers) {
1190  LogTrace("TrackExtenderWithMTD") << "Hit search: BTL layer at R= "
1191  << static_cast<const BarrelDetLayer*>(ilay)->specificSurface().radius();
1192 
1193  fillMatchingHits(ilay, tsos, traj, pmag2, pathlength0, trs0, hits, prop, bs, vtxTime, matchVertex, output, bestHit);
1194  }
1195 
1196  return output;
1197 }
1198 
1199 template <class TrackCollection>
1201  const TrajectoryStateOnSurface& tsos,
1202  const Trajectory& traj,
1203  const float pmag2,
1204  const float pathlength0,
1205  const TrackSegments& trs0,
1207  const MTDDetLayerGeometry* geo,
1208  const MagneticField* field,
1209  const Propagator* prop,
1210  const reco::BeamSpot& bs,
1211  const float vtxTime,
1212  const bool matchVertex,
1213  MTDHitMatchingInfo& bestHit) const {
1214  const vector<const DetLayer*>& layers = geo->allETLLayers();
1215 
1217  bestHit = MTDHitMatchingInfo();
1218  for (const DetLayer* ilay : layers) {
1219  const BoundDisk& disk = static_cast<const ForwardDetLayer*>(ilay)->specificSurface();
1220  const float diskZ = disk.position().z();
1221 
1222  if (tsos.globalPosition().z() * diskZ < 0)
1223  continue; // only propagate to the disk that's on the same side
1224 
1225  LogTrace("TrackExtenderWithMTD") << "Hit search: ETL disk at Z = " << diskZ;
1226 
1227  fillMatchingHits(ilay, tsos, traj, pmag2, pathlength0, trs0, hits, prop, bs, vtxTime, matchVertex, output, bestHit);
1228  }
1229 
1230  // the ETL hits order must be from the innermost to the outermost
1231 
1232  if (output.size() == 2) {
1233  if (std::abs(output[0]->globalPosition().z()) > std::abs(output[1]->globalPosition().z())) {
1234  std::reverse(output.begin(), output.end());
1235  }
1236  }
1237  return output;
1238 }
1239 
1240 template <class TrackCollection>
1242  const TrajectoryStateOnSurface& tsos,
1243  const Trajectory& traj,
1244  const float pmag2,
1245  const float pathlength0,
1246  const TrackSegments& trs0,
1248  const Propagator* prop,
1249  const reco::BeamSpot& bs,
1250  const float& vtxTime,
1251  const bool matchVertex,
1253  MTDHitMatchingInfo& bestHit) const {
1254  std::set<MTDHitMatchingInfo> hitsInLayer;
1255  bool hitMatched = false;
1256 
1257  using namespace std::placeholders;
1258  auto find_hits = std::bind(find_hits_in_dets,
1259  std::cref(hits),
1260  std::cref(traj),
1261  ilay,
1262  std::cref(tsos),
1263  pmag2,
1264  pathlength0,
1265  trs0,
1266  _1,
1267  std::cref(bs),
1268  bsTimeSpread_,
1269  prop,
1270  theEstimator.get(),
1271  _2,
1272  std::ref(hitsInLayer));
1273 
1274  if (useVertex_ && matchVertex)
1275  find_hits(vtxTime, true);
1276  else
1277  find_hits(0, false);
1278 
1279  float spaceChi2Cut = ilay->isBarrel() ? btlChi2Cut_ : etlChi2Cut_;
1280  float timeChi2Cut = ilay->isBarrel() ? btlTimeChi2Cut_ : etlTimeChi2Cut_;
1281 
1282  //just take the first hit because the hits are sorted on their matching quality
1283  if (!hitsInLayer.empty()) {
1284  //check hits to pass minimum quality matching requirements
1285  auto const& firstHit = *hitsInLayer.begin();
1286  LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: matching trial 1: estChi2= " << firstHit.estChi2
1287  << " timeChi2= " << firstHit.timeChi2;
1288  if (firstHit.estChi2 < spaceChi2Cut && firstHit.timeChi2 < timeChi2Cut) {
1289  hitMatched = true;
1290  output.push_back(hitbuilder_->build(firstHit.hit));
1291  if (firstHit < bestHit)
1292  bestHit = firstHit;
1293  }
1294  }
1295 
1296  if (useVertex_ && matchVertex && !hitMatched) {
1297  //try a second search with beamspot hypothesis
1298  hitsInLayer.clear();
1299  find_hits(0, false);
1300  if (!hitsInLayer.empty()) {
1301  auto const& firstHit = *hitsInLayer.begin();
1302  LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: matching trial 2: estChi2= " << firstHit.estChi2
1303  << " timeChi2= " << firstHit.timeChi2;
1304  if (firstHit.timeChi2 < timeChi2Cut) {
1305  if (firstHit.estChi2 < spaceChi2Cut) {
1306  hitMatched = true;
1307  output.push_back(hitbuilder_->build(firstHit.hit));
1308  if (firstHit < bestHit)
1309  bestHit = firstHit;
1310  }
1311  }
1312  }
1313  }
1314 
1315 #ifdef EDM_ML_DEBUG
1316  if (hitMatched) {
1317  LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: matched hit with time: " << bestHit.hit->time()
1318  << " +/- " << bestHit.hit->timeError();
1319  } else {
1320  LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: no matched hit";
1321  }
1322 #endif
1323 }
1324 
1325 //below is unfortunately ripped from other places but
1326 //since track producer doesn't know about MTD we have to do this
1327 template <class TrackCollection>
1329  const Trajectory& traj,
1330  const Trajectory& trajWithMtd,
1331  const reco::BeamSpot& bs,
1332  const MagneticField* field,
1333  const Propagator* thePropagator,
1334  bool hasMTD,
1335  float& pathLengthOut,
1336  float& tmtdOut,
1337  float& sigmatmtdOut,
1338  GlobalPoint& tmtdPosOut,
1339  float& tofpi,
1340  float& tofk,
1341  float& tofp,
1342  float& sigmatofpi,
1343  float& sigmatofk,
1344  float& sigmatofp) const {
1346  bool tsbcl_status = getTrajectoryStateClosestToBeamLine(traj, bs, thePropagator, tscbl);
1347 
1348  if (!tsbcl_status)
1349  return reco::Track();
1350 
1351  GlobalPoint v = tscbl.trackStateAtPCA().position();
1352  math::XYZPoint pos(v.x(), v.y(), v.z());
1353  GlobalVector p = tscbl.trackStateAtPCA().momentum();
1354  math::XYZVector mom(p.x(), p.y(), p.z());
1355 
1356  int ndof = traj.ndof();
1357 
1358  float t0 = 0.f;
1359  float covt0t0 = -1.f;
1360  pathLengthOut = -1.f; // if there is no MTD flag the pathlength with -1
1361  tmtdOut = 0.f;
1362  sigmatmtdOut = -1.f;
1363  float betaOut = 0.f;
1364  float covbetabeta = -1.f;
1365 
1366  auto routput = [&]() {
1367  return reco::Track(traj.chiSquared(),
1368  int(ndof),
1369  pos,
1370  mom,
1371  tscbl.trackStateAtPCA().charge(),
1373  orig->algo(),
1375  t0,
1376  betaOut,
1377  covt0t0,
1378  covbetabeta);
1379  };
1380 
1381  //compute path length for time backpropagation, using first MTD hit for the momentum
1382  if (hasMTD) {
1383  float pathlength;
1384  TrackSegments trs;
1385  bool validpropagation = trackPathLength(trajWithMtd, bs, thePropagator, pathlength, trs);
1386  float thit = 0.f;
1387  float thiterror = -1.f;
1388  GlobalPoint thitpos{0., 0., 0.};
1389  bool validmtd = false;
1390 
1391  if (!validpropagation) {
1392  return routput();
1393  }
1394 
1395  uint32_t ihitcount(0), ietlcount(0);
1396  for (auto const& hit : trajWithMtd.measurements()) {
1397  if (hit.recHit()->geographicalId().det() == DetId::Forward &&
1398  ForwardSubdetector(hit.recHit()->geographicalId().subdetId()) == FastTime) {
1399  ihitcount++;
1400  if (MTDDetId(hit.recHit()->geographicalId()).mtdSubDetector() == MTDDetId::MTDType::ETL) {
1401  ietlcount++;
1402  }
1403  }
1404  }
1405 
1406  LogTrace("TrackExtenderWithMTD") << "TrackExtenderWithMTD: selected #hits " << ihitcount << " from ETL "
1407  << ietlcount;
1408 
1409  auto ihit1 = trajWithMtd.measurements().cbegin();
1410  if (ihitcount == 1) {
1411  const MTDTrackingRecHit* mtdhit = static_cast<const MTDTrackingRecHit*>((*ihit1).recHit()->hit());
1412  thit = mtdhit->time();
1413  thiterror = mtdhit->timeError();
1414  thitpos = mtdhit->globalPosition();
1415  validmtd = true;
1416  } else if (ihitcount == 2 && ietlcount == 2) {
1417  std::pair<float, float> lastStep = trs.getSegmentPathAndMom2(0);
1418  float etlpathlength = std::abs(lastStep.first * c_cm_ns);
1419  //
1420  // The information of the two ETL hits is combined and attributed to the innermost hit
1421  //
1422  if (etlpathlength == 0.f) {
1423  validpropagation = false;
1424  } else {
1425  pathlength -= etlpathlength;
1426  trs.removeFirstSegment();
1427  const MTDTrackingRecHit* mtdhit1 = static_cast<const MTDTrackingRecHit*>((*ihit1).recHit()->hit());
1428  const MTDTrackingRecHit* mtdhit2 = static_cast<const MTDTrackingRecHit*>((*(ihit1 + 1)).recHit()->hit());
1429  TrackTofPidInfo tofInfo = computeTrackTofPidInfo(
1430  lastStep.second, etlpathlength, trs, mtdhit1->time(), mtdhit1->timeError(), 0.f, 0.f, true, TofCalc::kCost);
1431  //
1432  // Protect against incompatible times
1433  //
1434  float err1 = tofInfo.dterror2;
1435  float err2 = mtdhit2->timeError() * mtdhit2->timeError();
1436  if (cms_rounding::roundIfNear0(err1) == 0.f || cms_rounding::roundIfNear0(err2) == 0.f) {
1437  edm::LogError("TrackExtenderWithMTD")
1438  << "MTD tracking hits with zero time uncertainty: " << err1 << " " << err2;
1439  } else {
1440  if ((tofInfo.dt - mtdhit2->time()) * (tofInfo.dt - mtdhit2->time()) < (err1 + err2) * etlTimeChi2Cut_) {
1441  //
1442  // Subtract the ETL time of flight from the outermost measurement, and combine it in a weighted average with the innermost
1443  // the mass ambiguity related uncertainty on the time of flight is added as an additional uncertainty
1444  //
1445  err1 = 1.f / err1;
1446  err2 = 1.f / err2;
1447  thiterror = 1.f / (err1 + err2);
1448  thit = (tofInfo.dt * err1 + mtdhit2->time() * err2) * thiterror;
1449  thiterror = std::sqrt(thiterror);
1450  thitpos = mtdhit2->globalPosition();
1451  LogTrace("TrackExtenderWithMTD")
1452  << "TrackExtenderWithMTD: p trk = " << p.mag() << " ETL hits times/errors: 1) " << mtdhit1->time()
1453  << " +/- " << mtdhit1->timeError() << " , 2) " << mtdhit2->time() << " +/- " << mtdhit2->timeError()
1454  << " extrapolated time1: " << tofInfo.dt << " +/- " << tofInfo.dterror << " average = " << thit
1455  << " +/- " << thiterror << "\n hit1 pos: " << mtdhit1->globalPosition()
1456  << " hit2 pos: " << mtdhit2->globalPosition() << " etl path length " << etlpathlength << std::endl;
1457  validmtd = true;
1458  } else {
1459  // if back extrapolated time of the outermost measurement not compatible with the innermost, keep the one with smallest error
1460  if (err1 <= err2) {
1461  thit = tofInfo.dt;
1462  thiterror = tofInfo.dterror;
1463  validmtd = true;
1464  } else {
1465  thit = mtdhit2->time();
1466  thiterror = mtdhit2->timeError();
1467  validmtd = true;
1468  }
1469  }
1470  }
1471  }
1472  } else {
1473  edm::LogInfo("TrackExtenderWithMTD")
1474  << "MTD hits #" << ihitcount << "ETL hits #" << ietlcount << " anomalous pattern, skipping...";
1475  }
1476 
1477  if (validmtd && validpropagation) {
1478  //here add the PID uncertainty for later use in the 1st step of 4D vtx reconstruction
1479  TrackTofPidInfo tofInfo = computeTrackTofPidInfo(
1480  p.mag2(), pathlength, trs, thit, thiterror, 0.f, 0.f, true, TofCalc::kSegm, SigmaTofCalc::kCost);
1481 
1482  pathLengthOut = pathlength; // set path length if we've got a timing hit
1483  tmtdOut = thit;
1484  sigmatmtdOut = thiterror;
1485  tmtdPosOut = thitpos;
1486  t0 = tofInfo.dt;
1487  covt0t0 = tofInfo.dterror2;
1488  betaOut = tofInfo.beta_pi;
1489  covbetabeta = tofInfo.betaerror * tofInfo.betaerror;
1490  tofpi = tofInfo.dt_pi;
1491  tofk = tofInfo.dt_k;
1492  tofp = tofInfo.dt_p;
1493  sigmatofpi = tofInfo.sigma_dt_pi;
1494  sigmatofk = tofInfo.sigma_dt_k;
1495  sigmatofp = tofInfo.sigma_dt_p;
1496  }
1497  }
1498 
1499  return routput();
1500 }
1501 
1502 template <class TrackCollection>
1504  static const string metname = "TrackExtenderWithMTD";
1505 
1506  const Trajectory::RecHitContainer transRecHits = trajectory.recHits();
1507 
1508  // put the collection of TrackingRecHit in the event
1509 
1510  // sets the outermost and innermost TSOSs
1511  // ToDo: validation for track states with MTD
1512  TrajectoryStateOnSurface outerTSOS;
1513  TrajectoryStateOnSurface innerTSOS;
1514  unsigned int innerId = 0, outerId = 0;
1516  DetId outerDetId;
1517 
1518  if (trajectory.direction() == alongMomentum) {
1519  LogTrace(metname) << "alongMomentum";
1520  outerTSOS = trajectory.lastMeasurement().updatedState();
1521  innerTSOS = trajectory.firstMeasurement().updatedState();
1522  outerId = trajectory.lastMeasurement().recHit()->geographicalId().rawId();
1523  innerId = trajectory.firstMeasurement().recHit()->geographicalId().rawId();
1524  outerRecHit = trajectory.lastMeasurement().recHit();
1525  outerDetId = trajectory.lastMeasurement().recHit()->geographicalId();
1526  } else if (trajectory.direction() == oppositeToMomentum) {
1527  LogTrace(metname) << "oppositeToMomentum";
1528  outerTSOS = trajectory.firstMeasurement().updatedState();
1529  innerTSOS = trajectory.lastMeasurement().updatedState();
1530  outerId = trajectory.firstMeasurement().recHit()->geographicalId().rawId();
1531  innerId = trajectory.lastMeasurement().recHit()->geographicalId().rawId();
1532  outerRecHit = trajectory.firstMeasurement().recHit();
1533  outerDetId = trajectory.firstMeasurement().recHit()->geographicalId();
1534  } else
1535  LogError(metname) << "Wrong propagation direction!";
1536 
1537  const GeomDet* outerDet = gtg_->idToDet(outerDetId);
1538  GlobalPoint outerTSOSPos = outerTSOS.globalParameters().position();
1539  bool inside = outerDet->surface().bounds().inside(outerDet->toLocal(outerTSOSPos));
1540 
1541  GlobalPoint hitPos =
1542  (outerRecHit->isValid()) ? outerRecHit->globalPosition() : outerTSOS.globalParameters().position();
1543 
1544  if (!inside) {
1545  LogTrace(metname) << "The Global Muon outerMostMeasurementState is not compatible with the recHit detector!"
1546  << " Setting outerMost postition to recHit position if recHit isValid: "
1547  << outerRecHit->isValid();
1548  LogTrace(metname) << "From " << outerTSOSPos << " to " << hitPos;
1549  }
1550 
1551  //build the TrackExtra
1552  GlobalPoint v = (inside) ? outerTSOSPos : hitPos;
1553  GlobalVector p = outerTSOS.globalParameters().momentum();
1554  math::XYZPoint outpos(v.x(), v.y(), v.z());
1555  math::XYZVector outmom(p.x(), p.y(), p.z());
1556 
1557  v = innerTSOS.globalParameters().position();
1558  p = innerTSOS.globalParameters().momentum();
1559  math::XYZPoint inpos(v.x(), v.y(), v.z());
1560  math::XYZVector inmom(p.x(), p.y(), p.z());
1561 
1562  reco::TrackExtra trackExtra(outpos,
1563  outmom,
1564  true,
1565  inpos,
1566  inmom,
1567  true,
1568  outerTSOS.curvilinearError(),
1569  outerId,
1570  innerTSOS.curvilinearError(),
1571  innerId,
1572  trajectory.direction(),
1573  trajectory.seedRef());
1574 
1575  return trackExtra;
1576 }
1577 
1578 template <class TrackCollection>
1580  stringstream output;
1581 
1582  const BoundSurface* sur = nullptr;
1583  const BoundCylinder* bc = nullptr;
1584  const BoundDisk* bd = nullptr;
1585 
1586  sur = &(layer->surface());
1587  if ((bc = dynamic_cast<const BoundCylinder*>(sur))) {
1588  output << " Cylinder of radius: " << bc->radius() << endl;
1589  } else if ((bd = dynamic_cast<const BoundDisk*>(sur))) {
1590  output << " Disk at: " << bd->position().z() << endl;
1591  }
1592  return output.str();
1593 }
1594 
1595 //define this as a plug-in
1602 
size
Write out results.
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
float dt
Definition: AMPTWrapper.h:136
TrackExtenderWithMTDT< reco::TrackCollection > TrackExtenderWithMTD
edm::EDPutToken btlMatchTimeChi2Token_
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
bool isValid() const
Definition: Trajectory.h:257
edm::ESGetToken< TransientTrackingRecHitBuilder, TransientRecHitRecord > hitbuilderToken_
const CurvilinearTrajectoryError & curvilinearError() const
edm::EDPutToken btlMatchChi2Token_
std::vector< unsigned char > Chi2sFive
constexpr double c_cm_ns
virtual bool inside(const Local3DPoint &) const =0
Determine if the point is inside the bounds.
edm::EDPutToken tofpOrigTrkToken_
bool operator<(IOVSyncValue const &iLHS, IOVSyncValue const &iRHS)
Definition: IOVSyncValue.h:34
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > gtgToken_
edm::EDPutToken sigmatofkOrigTrkToken_
T z() const
Definition: PV3DBase.h:61
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
Definition: GlobalPoint.h:10
edm::EDPutToken etlMatchChi2Token_
edm::View< TrackType > InputCollection
edm::EDGetTokenT< GlobalPoint > genVtxPositionToken_
float chiSquared() const
Definition: Trajectory.h:241
edm::RefToBase< TrajectorySeed > seedRef(void) const
Definition: Trajectory.h:303
Log< level::Error, false > LogError
T mag2() const
Definition: PV3DBase.h:63
const SurfaceType & surface() const
edm::EDPutToken tofkOrigTrkToken_
edm::EDPutToken etlMatchTimeChi2Token_
ForwardSubdetector
edm::EDPutToken sigmat0OrigTrkToken_
void fillValueMap(edm::Event &iEvent, const H &handle, const std::vector< T > &vec, const edm::EDPutToken &token) const
edm::EDPutToken tmtdOrigTrkToken_
Definition: Electron.h:6
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:150
edm::EDGetTokenT< MTDTrackingDetSetVector > hitsToken_
#define LogTrace(id)
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
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
Definition: Trajectory.h:178
GlobalPoint position() const
edm::EDPutToken sigmatofpiOrigTrkToken_
string dumpLayer(const DetLayer *layer) const
std::vector< LocalTrajectoryParameters > TrajParams
int iEvent
Definition: GenABIO.cc:224
int ndof(bool bon=true) const
Definition: Trajectory.cc:97
GlobalPoint globalPosition() const
edm::EDPutToken tmtdPosOrigTrkToken_
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > builderToken_
const std::string mtdRecHitBuilder_
T sqrt(T t)
Definition: SSEVec.h:23
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
Definition: Trajectory.cc:133
TrackCharge charge() const
ConstRecHitContainer recHits() const
Definition: Trajectory.h:186
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, GlobalPoint &tmtdPosOut, float &tofpi, float &tofk, float &tofp, float &sigmatofpi, float &sigmatofk, float &sigmatofp) const
GlobalVector momentum() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDPutToken assocOrigTrkToken_
const std::string propagator_
edm::EDPutToken outermostHitPositionToken_
edm::ESGetToken< Propagator, TrackingComponentsRecord > propToken_
virtual std::pair< TrajectoryStateOnSurface, double > propagateWithPath(const FreeTrajectoryState &, const Surface &) const final
Definition: Propagator.cc:10
constexpr valType roundIfNear0(valType value, double tolerance=1.e-7)
Definition: Rounding.h:11
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
edm::ESHandle< TransientTrackingRecHitBuilder > hitbuilder_
const std::vector< const DetLayer * > & allBTLLayers() const
return the BTL DetLayers (barrel), inside-out
std::vector< ConstRecHitPointer > ConstRecHitContainer
constexpr double c_inv
const std::string transientTrackBuilder_
Log< level::Info, false > LogInfo
Definition: DetId.h:17
GlobalPoint globalPosition() const final
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:10
ConstRecHitContainer RecHitContainer
Definition: Trajectory.h:42
constexpr NumType convertMmToCm(NumType millimeters)
Definition: angle_units.h:44
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const CurvilinearTrajectoryError & curvilinearError() const
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
reco::TrackExtra buildTrackExtra(const Trajectory &trajectory) const
A 2D TrackerRecHit with time and time error information.
CurvilinearTrajectoryError curvilinearError(const PerigeeTrajectoryError &perigeeError, const GlobalTrajectoryParameters &gtp)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
TrackingRecHit::ConstRecHitPointer ConstRecHitPointer
edm::ESHandle< TransientTrackBuilder > builder_
float timeError() const
static void fillPSetDescription(edm::ParameterSetDescription &descriptions, bool doPredictionsOnly=false, const std::string &fitter="KFFitterForRefitInsideOut", const std::string &smoother="KFSmootherForRefitInsideOut", const std::string &propagator="SmartPropagatorAnyRK", const std::string &refitDirection="alongMomentum", bool refitRPCHits=true, const std::string &trackerRecHitBuilder="WithTrackAngle", const std::string &muonRecHitBuilder="MuonRecHitBuilder", const std::string &mtdRecHitBuilder="MTDRecHitBuilder")
fillDescriptions
std::unique_ptr< TrackTransformer > theTransformer
TrajectoryStateOnSurface const & updatedState() const
edm::EDPutToken npixBarrelToken_
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magfldToken_
edm::EDPutToken npixEndcapToken_
fixed size matrix
HLT enums.
TrajectoryMeasurement const & firstMeasurement() const
Definition: Trajectory.h:166
static int position[264][3]
Definition: ReadPGInfo.cc:289
TrajectoryMeasurement const & closestMeasurement(GlobalPoint) const
Definition: Trajectory.cc:256
edm::EDPutToken betaOrigTrkToken_
FreeTrajectoryState const * freeState(bool withErrors=true) const
ESTransientHandle< T > getTransientHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:141
edm::EDGetTokenT< TrajTrackAssociationCollection > trajTrackAToken_
float time() const
edm::EDPutToken pathLengthOrigTrkToken_
Definition: output.py:1
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_
ordering
Definition: config.py:7
edm::EDGetTokenT< VertexCollection > vtxToken_
if(threadIdxLocalY==0 &&threadIdxLocalX==0)
void produce(edm::Event &ev, const edm::EventSetup &es) final
bool isBarrel() const
Definition: DetLayer.h:31
TrackCollection::value_type TrackType
edm::EDGetTokenT< InputCollection > tracksToken_
const std::vector< const DetLayer * > & allETLLayers() const
return the ETL DetLayers (endcap), -Z to +Z
def move(src, dest)
Definition: eostools.py:511
RefitDirection::GeometricalDirection checkRecHitsOrdering(TransientTrackingRecHit::ConstRecHitContainer const &recHits) const
ConstRecHitPointer const & recHit() const
edm::EDPutToken sigmatofpOrigTrkToken_
#define LogDebug(id)
edm::ESGetToken< MTDDetLayerGeometry, MTDRecoGeometryRecord > dlgeoToken_
#define m_pi
Definition: RPCConst.cc:8
const Bounds & bounds() const
Definition: Surface.h:87