CMS 3D CMS Logo

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