CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TrackExtenderWithMTD.cc
Go to the documentation of this file.
8 
11 
14 
16 
18 
21 
26 
28 
32 
34 
37 
40 
42 
45 
46 #include <sstream>
47 
49 
52 
55 #include "CLHEP/Units/GlobalPhysicalConstants.h"
57 
60 
61 using namespace std;
62 using namespace edm;
63 using namespace reco;
64 
65 namespace {
66  class MTDHitMatchingInfo {
67  public:
68  MTDHitMatchingInfo() {
69  hit = nullptr;
72  }
73 
74  //Operator used to sort the hits while performing the matching step at the MTD
75  inline bool operator<(const MTDHitMatchingInfo& m2) const {
76  //only for good matching in time use estChi2, otherwise use mostly time compatibility
77  constexpr double chi2_cut = 10.;
78  constexpr double low_weight = 3.;
79  constexpr double high_weight = 8.;
80  if (timeChi2 < chi2_cut && m2.timeChi2 < chi2_cut)
81  return chi2(low_weight) < m2.chi2(low_weight);
82  else
83  return chi2(high_weight) < m2.chi2(high_weight);
84  }
85 
86  inline double chi2(double timeWeight = 1.) const { return estChi2 + timeWeight * timeChi2; }
87 
88  const MTDTrackingRecHit* hit;
89  double estChi2;
90  double timeChi2;
91  };
92 
93  struct TrackTofPidInfo {
94  double tmtd;
95  double tmtderror;
96  double pathlength;
97 
98  double betaerror;
99 
100  double dt;
101  double dterror;
102  double dtchi2;
103 
104  double dt_best;
105  double dterror_best;
106  double dtchi2_best;
107 
108  double gammasq_pi;
109  double beta_pi;
110  double dt_pi;
111 
112  double gammasq_k;
113  double beta_k;
114  double dt_k;
115 
116  double gammasq_p;
117  double beta_p;
118  double dt_p;
119 
120  double prob_pi;
121  double prob_k;
122  double prob_p;
123  };
124 
125  const TrackTofPidInfo computeTrackTofPidInfo(double magp2,
126  double length,
127  double t_mtd,
128  double t_mtderr,
129  double t_vtx,
130  double t_vtx_err,
131  bool addPIDError = true) {
132  constexpr double m_pi = 0.13957018;
133  constexpr double m_pi_inv2 = 1.0 / m_pi / m_pi;
134  constexpr double m_k = 0.493677;
135  constexpr double m_k_inv2 = 1.0 / m_k / m_k;
136  constexpr double m_p = 0.9382720813;
137  constexpr double m_p_inv2 = 1.0 / m_p / m_p;
138  constexpr double c_cm_ns = geant_units::operators::convertMmToCm(CLHEP::c_light); // [mm/ns] -> [cm/ns]
139  constexpr double c_inv = 1.0 / c_cm_ns;
140 
141  TrackTofPidInfo tofpid;
142 
143  tofpid.tmtd = t_mtd;
144  tofpid.tmtderror = t_mtderr;
145  tofpid.pathlength = length;
146 
147  tofpid.gammasq_pi = 1. + magp2 * m_pi_inv2;
148  tofpid.beta_pi = std::sqrt(1. - 1. / tofpid.gammasq_pi);
149  tofpid.dt_pi = tofpid.pathlength / tofpid.beta_pi * c_inv;
150 
151  tofpid.gammasq_k = 1. + magp2 * m_k_inv2;
152  tofpid.beta_k = std::sqrt(1. - 1. / tofpid.gammasq_k);
153  tofpid.dt_k = tofpid.pathlength / tofpid.beta_k * c_inv;
154 
155  tofpid.gammasq_p = 1. + magp2 * m_p_inv2;
156  tofpid.beta_p = std::sqrt(1. - 1. / tofpid.gammasq_p);
157  tofpid.dt_p = tofpid.pathlength / tofpid.beta_p * c_inv;
158 
159  tofpid.dt = tofpid.tmtd - tofpid.dt_pi - t_vtx; //assume by default the pi hypothesis
160  tofpid.dterror = sqrt(tofpid.tmtderror * tofpid.tmtderror + t_vtx_err * t_vtx_err);
161  tofpid.betaerror = 0;
162  if (addPIDError) {
163  tofpid.dterror =
164  sqrt(tofpid.dterror * tofpid.dterror + (tofpid.dt_p - tofpid.dt_pi) * (tofpid.dt_p - tofpid.dt_pi));
165  tofpid.betaerror = tofpid.beta_p - tofpid.beta_pi;
166  }
167 
168  tofpid.dtchi2 = (tofpid.dt * tofpid.dt) / (tofpid.dterror * tofpid.dterror);
169 
170  tofpid.dt_best = tofpid.dt;
171  tofpid.dterror_best = tofpid.dterror;
172  tofpid.dtchi2_best = tofpid.dtchi2;
173 
174  tofpid.prob_pi = -1.;
175  tofpid.prob_k = -1.;
176  tofpid.prob_p = -1.;
177 
178  if (!addPIDError) {
179  //*TODO* deal with heavier nucleons and/or BSM case here?
180  double chi2_pi = tofpid.dtchi2;
181  double chi2_k =
182  (tofpid.tmtd - tofpid.dt_k - t_vtx) * (tofpid.tmtd - tofpid.dt_k - t_vtx) / (tofpid.dterror * tofpid.dterror);
183  double chi2_p =
184  (tofpid.tmtd - tofpid.dt_p - t_vtx) * (tofpid.tmtd - tofpid.dt_p - t_vtx) / (tofpid.dterror * tofpid.dterror);
185 
186  double rawprob_pi = exp(-0.5 * chi2_pi);
187  double rawprob_k = exp(-0.5 * chi2_k);
188  double rawprob_p = exp(-0.5 * chi2_p);
189  double normprob = 1. / (rawprob_pi + rawprob_k + rawprob_p);
190 
191  tofpid.prob_pi = rawprob_pi * normprob;
192  tofpid.prob_k = rawprob_k * normprob;
193  tofpid.prob_p = rawprob_p * normprob;
194 
195  double prob_heavy = 1. - tofpid.prob_pi;
196  constexpr double heavy_threshold = 0.75;
197 
198  if (prob_heavy > heavy_threshold) {
199  if (chi2_k < chi2_p) {
200  tofpid.dt_best = (tofpid.tmtd - tofpid.dt_k - t_vtx);
201  tofpid.dtchi2_best = chi2_k;
202  } else {
203  tofpid.dt_best = (tofpid.tmtd - tofpid.dt_p - t_vtx);
204  tofpid.dtchi2_best = chi2_p;
205  }
206  }
207  }
208  return tofpid;
209  }
210 
211  bool getTrajectoryStateClosestToBeamLine(const Trajectory& traj,
212  const reco::BeamSpot& bs,
213  const Propagator* thePropagator,
215  // get the state closest to the beamline
216  TrajectoryStateOnSurface stateForProjectionToBeamLineOnSurface =
217  traj.closestMeasurement(GlobalPoint(bs.x0(), bs.y0(), bs.z0())).updatedState();
218 
219  if (!stateForProjectionToBeamLineOnSurface.isValid()) {
220  edm::LogError("CannotPropagateToBeamLine") << "the state on the closest measurement isnot valid. skipping track.";
221  return false;
222  }
223 
224  const FreeTrajectoryState& stateForProjectionToBeamLine = *stateForProjectionToBeamLineOnSurface.freeState();
225 
226  TSCBLBuilderWithPropagator tscblBuilder(*thePropagator);
227  tscbl = tscblBuilder(stateForProjectionToBeamLine, bs);
228 
229  return tscbl.isValid();
230  }
231 
232  bool trackPathLength(const Trajectory& traj,
234  const Propagator* thePropagator,
235  double& pathlength) {
236  pathlength = 0.;
237 
238  bool validpropagation = true;
239  double pathlength1 = 0.;
240  double pathlength2 = 0.;
241 
242  //add pathlength layer by layer
243  for (auto it = traj.measurements().begin(); it != traj.measurements().end() - 1; ++it) {
244  const auto& propresult = thePropagator->propagateWithPath(it->updatedState(), (it + 1)->updatedState().surface());
245  double layerpathlength = std::abs(propresult.second);
246  if (layerpathlength == 0.) {
247  validpropagation = false;
248  }
249  pathlength1 += layerpathlength;
250  }
251 
252  //add distance from bs to first measurement
253  auto const& tscblPCA = tscbl.trackStateAtPCA();
254  auto const& aSurface = traj.direction() == alongMomentum ? traj.firstMeasurement().updatedState().surface()
255  : traj.lastMeasurement().updatedState().surface();
256  pathlength2 = thePropagator->propagateWithPath(tscblPCA, aSurface).second;
257  if (pathlength2 == 0.) {
258  validpropagation = false;
259  }
260  pathlength = pathlength1 + pathlength2;
261 
262  return validpropagation;
263  }
264 
265  bool trackPathLength(const Trajectory& traj,
266  const reco::BeamSpot& bs,
267  const Propagator* thePropagator,
268  double& pathlength) {
269  pathlength = 0.;
270 
272  bool tscbl_status = getTrajectoryStateClosestToBeamLine(traj, bs, thePropagator, tscbl);
273 
274  if (!tscbl_status)
275  return false;
276 
277  return trackPathLength(traj, tscbl, thePropagator, pathlength);
278  }
279 
280 } // namespace
281 
282 template <class TrackCollection>
284 public:
287 
289 
290  template <class H, class T>
291  void fillValueMap(edm::Event& iEvent, const H& handle, const std::vector<T>& vec, const edm::EDPutToken& token) const;
292 
293  void produce(edm::Event& ev, const edm::EventSetup& es) final;
294 
295  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
296 
298  const Trajectory& traj,
299  const double,
300  const double,
302  const MTDDetLayerGeometry*,
303  const MagneticField* field,
304  const Propagator* prop,
305  const reco::BeamSpot& bs,
306  const double vtxTime,
307  const bool matchVertex,
308  MTDHitMatchingInfo& bestHit) const;
309 
311  const Trajectory& traj,
312  const double,
313  const double,
315  const MTDDetLayerGeometry*,
316  const MagneticField* field,
317  const Propagator* prop,
318  const reco::BeamSpot& bs,
319  const double vtxTime,
320  const bool matchVertex,
321  MTDHitMatchingInfo& bestHit) const;
322 
323  void fillMatchingHits(const DetLayer*,
325  const Trajectory&,
326  const double,
327  const double,
329  const Propagator*,
330  const reco::BeamSpot&,
331  const double&,
332  const bool,
334  MTDHitMatchingInfo&) const;
335 
338  if (!recHits.empty()) {
339  GlobalPoint first = gtg_->idToDet(recHits.front()->geographicalId())->position();
340  GlobalPoint last = gtg_->idToDet(recHits.back()->geographicalId())->position();
341 
342  // maybe perp2?
343  auto rFirst = first.mag2();
344  auto rLast = last.mag2();
345  if (rFirst < rLast)
347  if (rFirst > rLast)
349  }
350  LogDebug("TrackExtenderWithMTD") << "Impossible to determine the rechits order" << endl;
352  }
353 
354  reco::Track buildTrack(const reco::Track&,
355  const Trajectory&,
356  const Trajectory&,
357  const reco::BeamSpot&,
358  const MagneticField* field,
359  const Propagator* prop,
360  bool hasMTD,
361  float& pathLength,
362  float& tmtdOut,
363  float& sigmatmtdOut) const;
364  reco::TrackExtra buildTrackExtra(const Trajectory& trajectory) const;
365 
366  string dumpLayer(const DetLayer* layer) const;
367 
368 private:
383 
390 
391  const bool updateTraj_, updateExtra_, updatePattern_;
392  const std::string mtdRecHitBuilder_, propagator_, transientTrackBuilder_;
393  std::unique_ptr<MeasurementEstimator> theEstimator;
394  std::unique_ptr<TrackTransformer> theTransformer;
401 
406 
407  const float estMaxChi2_;
408  const float estMaxNSigma_;
409  const float btlChi2Cut_;
410  const float btlTimeChi2Cut_;
411  const float etlChi2Cut_;
412  const float etlTimeChi2Cut_;
413 
414  const bool useVertex_;
415  const bool useSimVertex_;
416  const float dzCut_;
417  const float bsTimeSpread_;
418 };
419 
420 template <class TrackCollection>
422  : tracksToken_(consumes<InputCollection>(iConfig.getParameter<edm::InputTag>("tracksSrc"))),
423  hitsToken_(consumes<MTDTrackingDetSetVector>(iConfig.getParameter<edm::InputTag>("hitsSrc"))),
424  bsToken_(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpotSrc"))),
425  updateTraj_(iConfig.getParameter<bool>("updateTrackTrajectory")),
426  updateExtra_(iConfig.getParameter<bool>("updateTrackExtra")),
427  updatePattern_(iConfig.getParameter<bool>("updateTrackHitPattern")),
428  mtdRecHitBuilder_(iConfig.getParameter<std::string>("MTDRecHitBuilder")),
429  propagator_(iConfig.getParameter<std::string>("Propagator")),
430  transientTrackBuilder_(iConfig.getParameter<std::string>("TransientTrackBuilder")),
431  estMaxChi2_(iConfig.getParameter<double>("estimatorMaxChi2")),
432  estMaxNSigma_(iConfig.getParameter<double>("estimatorMaxNSigma")),
433  btlChi2Cut_(iConfig.getParameter<double>("btlChi2Cut")),
434  btlTimeChi2Cut_(iConfig.getParameter<double>("btlTimeChi2Cut")),
435  etlChi2Cut_(iConfig.getParameter<double>("etlChi2Cut")),
436  etlTimeChi2Cut_(iConfig.getParameter<double>("etlTimeChi2Cut")),
437  useVertex_(iConfig.getParameter<bool>("useVertex")),
438  useSimVertex_(iConfig.getParameter<bool>("useSimVertex")),
439  dzCut_(iConfig.getParameter<double>("dZCut")),
440  bsTimeSpread_(iConfig.getParameter<double>("bsTimeSpread")) {
441  if (useVertex_) {
442  if (useSimVertex_) {
443  genVtxPositionToken_ = consumes<GlobalPoint>(iConfig.getParameter<edm::InputTag>("genVtxPositionSrc"));
444  genVtxTimeToken_ = consumes<float>(iConfig.getParameter<edm::InputTag>("genVtxTimeSrc"));
445  } else
446  vtxToken_ = consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("vtxSrc"));
447  }
448 
449  theEstimator = std::make_unique<Chi2MeasurementEstimator>(estMaxChi2_, estMaxNSigma_);
450  theTransformer = std::make_unique<TrackTransformer>(iConfig.getParameterSet("TrackTransformer"), consumesCollector());
451 
452  btlMatchChi2Token = produces<edm::ValueMap<float>>("btlMatchChi2");
453  etlMatchChi2Token = produces<edm::ValueMap<float>>("etlMatchChi2");
454  btlMatchTimeChi2Token = produces<edm::ValueMap<float>>("btlMatchTimeChi2");
455  etlMatchTimeChi2Token = produces<edm::ValueMap<float>>("etlMatchTimeChi2");
456  npixBarrelToken = produces<edm::ValueMap<int>>("npixBarrel");
457  npixEndcapToken = produces<edm::ValueMap<int>>("npixEndcap");
458  pOrigTrkToken = produces<edm::ValueMap<float>>("generalTrackp");
459  betaOrigTrkToken = produces<edm::ValueMap<float>>("generalTrackBeta");
460  t0OrigTrkToken = produces<edm::ValueMap<float>>("generalTrackt0");
461  sigmat0OrigTrkToken = produces<edm::ValueMap<float>>("generalTracksigmat0");
462  pathLengthOrigTrkToken = produces<edm::ValueMap<float>>("generalTrackPathLength");
463  tmtdOrigTrkToken = produces<edm::ValueMap<float>>("generalTracktmtd");
464  sigmatmtdOrigTrkToken = produces<edm::ValueMap<float>>("generalTracksigmatmtd");
465  assocOrigTrkToken = produces<edm::ValueMap<int>>("generalTrackassoc");
466 
467  builderToken_ = esConsumes<TransientTrackBuilder, TransientTrackRecord>(edm::ESInputTag("", transientTrackBuilder_));
469  esConsumes<TransientTrackingRecHitBuilder, TransientRecHitRecord>(edm::ESInputTag("", mtdRecHitBuilder_));
470  gtgToken_ = esConsumes<GlobalTrackingGeometry, GlobalTrackingGeometryRecord>();
471  dlgeoToken_ = esConsumes<MTDDetLayerGeometry, MTDRecoGeometryRecord>();
472  magfldToken_ = esConsumes<MagneticField, IdealMagneticFieldRecord>();
473  propToken_ = esConsumes<Propagator, TrackingComponentsRecord>(edm::ESInputTag("", propagator_));
474  ttopoToken_ = esConsumes<TrackerTopology, TrackerTopologyRcd>();
475 
476  produces<edm::OwnVector<TrackingRecHit>>();
477  produces<reco::TrackExtraCollection>();
478  produces<TrackCollection>();
479 }
480 
481 template <class TrackCollection>
484  desc.add<edm::InputTag>("tracksSrc", edm::InputTag("generalTracks"));
485  desc.add<edm::InputTag>("hitsSrc", edm::InputTag("mtdTrackingRecHits"));
486  desc.add<edm::InputTag>("beamSpotSrc", edm::InputTag("offlineBeamSpot"));
487  desc.add<edm::InputTag>("genVtxPositionSrc", edm::InputTag("genParticles:xyz0"));
488  desc.add<edm::InputTag>("genVtxTimeSrc", edm::InputTag("genParticles:t0"));
489  desc.add<edm::InputTag>("vtxSrc", edm::InputTag("offlinePrimaryVertices4D"));
490  desc.add<bool>("updateTrackTrajectory", true);
491  desc.add<bool>("updateTrackExtra", true);
492  desc.add<bool>("updateTrackHitPattern", true);
493  desc.add<std::string>("TransientTrackBuilder", "TransientTrackBuilder");
494  desc.add<std::string>("MTDRecHitBuilder", "MTDRecHitBuilder");
495  desc.add<std::string>("Propagator", "PropagatorWithMaterialForMTD");
497  false,
498  "KFFitterForRefitInsideOut",
499  "KFSmootherForRefitInsideOut",
500  "PropagatorWithMaterialForMTD",
501  "alongMomentum",
502  true,
503  "WithTrackAngle",
504  "MuonRecHitBuilder",
505  "MTDRecHitBuilder");
506  desc.add<edm::ParameterSetDescription>("TrackTransformer", transDesc);
507  desc.add<double>("estimatorMaxChi2", 500.);
508  desc.add<double>("estimatorMaxNSigma", 10.);
509  desc.add<double>("btlChi2Cut", 50.);
510  desc.add<double>("btlTimeChi2Cut", 10.);
511  desc.add<double>("etlChi2Cut", 50.);
512  desc.add<double>("etlTimeChi2Cut", 10.);
513  desc.add<bool>("useVertex", false);
514  desc.add<bool>("useSimVertex", false);
515  desc.add<double>("dZCut", 0.1);
516  desc.add<double>("bsTimeSpread", 0.2);
517  descriptions.add("trackExtenderWithMTDBase", desc);
518 }
519 
520 template <class TrackCollection>
521 template <class H, class T>
523  const H& handle,
524  const std::vector<T>& vec,
525  const edm::EDPutToken& token) const {
526  auto out = std::make_unique<edm::ValueMap<T>>();
527  typename edm::ValueMap<T>::Filler filler(*out);
528  filler.insert(handle, vec.begin(), vec.end());
529  filler.fill();
530  iEvent.put(token, std::move(out));
531 }
532 
533 template <class TrackCollection>
535  //this produces pieces of the track extra
536  Traj2TrackHits t2t;
537 
538  theTransformer->setServices(es);
539 
542 
543  gtg_ = es.getHandle(gtgToken_);
544 
545  auto geo = es.getTransientHandle(dlgeoToken_);
546 
547  auto magfield = es.getTransientHandle(magfldToken_);
548 
549  builder_ = es.getHandle(builderToken_);
550  hitbuilder_ = es.getHandle(hitbuilderToken_);
551 
552  auto propH = es.getTransientHandle(propToken_);
553  const Propagator* prop = propH.product();
554 
555  auto httopo = es.getTransientHandle(ttopoToken_);
556  const TrackerTopology& ttopo = *httopo;
557 
558  auto output = std::make_unique<TrackCollection>();
559  auto extras = std::make_unique<reco::TrackExtraCollection>();
560  auto outhits = std::make_unique<edm::OwnVector<TrackingRecHit>>();
561 
562  std::vector<float> btlMatchChi2;
563  std::vector<float> etlMatchChi2;
564  std::vector<float> btlMatchTimeChi2;
565  std::vector<float> etlMatchTimeChi2;
566  std::vector<int> npixBarrel;
567  std::vector<int> npixEndcap;
568  std::vector<float> pOrigTrkRaw;
569  std::vector<float> betaOrigTrkRaw;
570  std::vector<float> t0OrigTrkRaw;
571  std::vector<float> sigmat0OrigTrkRaw;
572  std::vector<float> pathLengthsOrigTrkRaw;
573  std::vector<float> tmtdOrigTrkRaw;
574  std::vector<float> sigmatmtdOrigTrkRaw;
575  std::vector<int> assocOrigTrkRaw;
576 
577  auto const tracksH = ev.getHandle(tracksToken_);
578  const auto& tracks = *tracksH;
579 
580  //MTD hits DetSet
581  const auto& hits = ev.get(hitsToken_);
582 
583  //beam spot
584  const auto& bs = ev.get(bsToken_);
585 
586  const Vertex* pv = nullptr;
587  if (useVertex_ && !useSimVertex_) {
588  auto const& vtxs = ev.get(vtxToken_);
589  if (!vtxs.empty())
590  pv = &vtxs[0];
591  }
592 
593  std::unique_ptr<math::XYZTLorentzVectorF> genPV(nullptr);
594  if (useVertex_ && useSimVertex_) {
595  const auto& genVtxPosition = ev.get(genVtxPositionToken_);
596  const auto& genVtxTime = ev.get(genVtxTimeToken_);
597  genPV = std::make_unique<math::XYZTLorentzVectorF>(
598  genVtxPosition.x(), genVtxPosition.y(), genVtxPosition.z(), genVtxTime);
599  }
600 
601  double vtxTime = 0.;
602  if (useVertex_) {
603  if (useSimVertex_ && genPV) {
604  vtxTime = genPV->t();
605  } else if (pv)
606  vtxTime = pv->t(); //already in ns
607  }
608 
609  std::vector<unsigned> track_indices;
610  unsigned itrack = 0;
611 
612  for (const auto& track : tracks) {
613  double trackVtxTime = 0.;
614  if (useVertex_) {
615  double dz;
616  if (useSimVertex_)
617  dz = std::abs(track.dz(math::XYZPoint(*genPV)));
618  else
619  dz = std::abs(track.dz(pv->position()));
620 
621  if (dz < dzCut_)
622  trackVtxTime = vtxTime;
623  }
624 
625  reco::TransientTrack ttrack(track, magfield.product(), gtg_);
626  const auto& trajs = theTransformer->transform(track);
627  auto thits = theTransformer->getTransientRecHits(ttrack);
629  MTDHitMatchingInfo mBTL, mETL;
630  if (!trajs.empty()) {
631  // get the outermost trajectory point on the track
632  TrajectoryStateOnSurface tsos = builder_->build(track).outermostMeasurementState();
634  bool tscbl_status = getTrajectoryStateClosestToBeamLine(trajs.front(), bs, prop, tscbl);
635 
636  if (tscbl_status) {
637  double pmag2 = tscbl.trackStateAtPCA().momentum().mag2();
638  double pathlength0;
639  trackPathLength(trajs.front(), tscbl, prop, pathlength0);
640 
641  const auto& btlhits = tryBTLLayers(tsos,
642  trajs.front(),
643  pmag2,
644  pathlength0,
645  hits,
646  geo.product(),
647  magfield.product(),
648  prop,
649  bs,
650  trackVtxTime,
651  trackVtxTime != 0.,
652  mBTL);
653  mtdthits.insert(mtdthits.end(), btlhits.begin(), btlhits.end());
654 
655  // in the future this should include an intermediate refit before propagating to the ETL
656  // for now it is ok
657  const auto& etlhits = tryETLLayers(tsos,
658  trajs.front(),
659  pmag2,
660  pathlength0,
661  hits,
662  geo.product(),
663  magfield.product(),
664  prop,
665  bs,
666  trackVtxTime,
667  trackVtxTime != 0.,
668  mETL);
669  mtdthits.insert(mtdthits.end(), etlhits.begin(), etlhits.end());
670  }
671  }
672 
673  auto ordering = checkRecHitsOrdering(thits);
675  thits.insert(thits.end(), mtdthits.begin(), mtdthits.end());
676  } else {
677  std::reverse(mtdthits.begin(), mtdthits.end());
678  mtdthits.insert(mtdthits.end(), thits.begin(), thits.end());
679  thits.swap(mtdthits);
680  }
681 
682  const auto& trajwithmtd = mtdthits.empty() ? trajs : theTransformer->transform(ttrack, thits);
683  float pMap = 0.f, betaMap = 0.f, t0Map = 0.f, sigmat0Map = -1.f, pathLengthMap = -1.f, tmtdMap = 0.f,
684  sigmatmtdMap = -1.f;
685  int iMap = -1;
686 
687  for (const auto& trj : trajwithmtd) {
688  const auto& thetrj = (updateTraj_ ? trj : trajs.front());
689  float pathLength = 0.f, tmtd = 0.f, sigmatmtd = -1.f;
690  reco::Track result = buildTrack(track,
691  thetrj,
692  trj,
693  bs,
694  magfield.product(),
695  prop,
696  !trajwithmtd.empty() && !mtdthits.empty(),
697  pathLength,
698  tmtd,
699  sigmatmtd);
700  if (result.ndof() >= 0) {
702  reco::TrackExtra::TrajParams trajParams;
704  size_t hitsstart = outhits->size();
705  if (updatePattern_) {
706  t2t(trj, *outhits, trajParams, chi2s); // this fills the output hit collection
707  } else {
708  t2t(thetrj, *outhits, trajParams, chi2s);
709  }
710  size_t hitsend = outhits->size();
711  extras->push_back(buildTrackExtra(trj)); // always push back the fully built extra, update by setting in track
712  extras->back().setHits(hitsRefProd, hitsstart, hitsend - hitsstart);
713  extras->back().setTrajParams(trajParams, chi2s);
714  //create the track
715  output->push_back(result);
716  btlMatchChi2.push_back(mBTL.hit ? mBTL.estChi2 : -1);
717  etlMatchChi2.push_back(mETL.hit ? mETL.estChi2 : -1);
718  btlMatchTimeChi2.push_back(mBTL.hit ? mBTL.timeChi2 : -1);
719  etlMatchTimeChi2.push_back(mETL.hit ? mETL.timeChi2 : -1);
720  pathLengthMap = pathLength;
721  tmtdMap = tmtd;
722  sigmatmtdMap = sigmatmtd;
723  auto& backtrack = output->back();
724  iMap = output->size() - 1;
725  pMap = backtrack.p();
726  betaMap = backtrack.beta();
727  t0Map = backtrack.t0();
728  sigmat0Map = std::copysign(std::sqrt(std::abs(backtrack.covt0t0())), backtrack.covt0t0());
729  reco::TrackExtraRef extraRef(extrasRefProd, extras->size() - 1);
730  backtrack.setExtra((updateExtra_ ? extraRef : track.extra()));
731  for (unsigned ihit = hitsstart; ihit < hitsend; ++ihit) {
732  backtrack.appendHitPattern((*outhits)[ihit], ttopo);
733  }
734  npixBarrel.push_back(backtrack.hitPattern().numberOfValidPixelBarrelHits());
735  npixEndcap.push_back(backtrack.hitPattern().numberOfValidPixelEndcapHits());
736  } else {
737  LogTrace("TrackExtenderWithMTD") << "Error in the MTD track refitting. This should not happen";
738  }
739  }
740 
741  pOrigTrkRaw.push_back(pMap);
742  betaOrigTrkRaw.push_back(betaMap);
743  t0OrigTrkRaw.push_back(t0Map);
744  sigmat0OrigTrkRaw.push_back(sigmat0Map);
745  pathLengthsOrigTrkRaw.push_back(pathLengthMap);
746  tmtdOrigTrkRaw.push_back(tmtdMap);
747  sigmatmtdOrigTrkRaw.push_back(sigmatmtdMap);
748  assocOrigTrkRaw.push_back(iMap);
749 
750  if (iMap == -1) {
751  btlMatchChi2.push_back(-1.);
752  etlMatchChi2.push_back(-1.);
753  btlMatchTimeChi2.push_back(-1.);
754  etlMatchTimeChi2.push_back(-1.);
755  npixBarrel.push_back(-1.);
756  npixEndcap.push_back(-1.);
757  }
758 
759  ++itrack;
760  }
761 
762  auto outTrksHandle = ev.put(std::move(output));
763  ev.put(std::move(extras));
764  ev.put(std::move(outhits));
765 
766  fillValueMap(ev, tracksH, btlMatchChi2, btlMatchChi2Token);
767  fillValueMap(ev, tracksH, etlMatchChi2, etlMatchChi2Token);
768  fillValueMap(ev, tracksH, btlMatchTimeChi2, btlMatchTimeChi2Token);
769  fillValueMap(ev, tracksH, etlMatchTimeChi2, etlMatchTimeChi2Token);
770  fillValueMap(ev, tracksH, npixBarrel, npixBarrelToken);
771  fillValueMap(ev, tracksH, npixEndcap, npixEndcapToken);
772  fillValueMap(ev, tracksH, pOrigTrkRaw, pOrigTrkToken);
773  fillValueMap(ev, tracksH, betaOrigTrkRaw, betaOrigTrkToken);
774  fillValueMap(ev, tracksH, t0OrigTrkRaw, t0OrigTrkToken);
775  fillValueMap(ev, tracksH, sigmat0OrigTrkRaw, sigmat0OrigTrkToken);
776  fillValueMap(ev, tracksH, pathLengthsOrigTrkRaw, pathLengthOrigTrkToken);
777  fillValueMap(ev, tracksH, tmtdOrigTrkRaw, tmtdOrigTrkToken);
778  fillValueMap(ev, tracksH, sigmatmtdOrigTrkRaw, sigmatmtdOrigTrkToken);
779  fillValueMap(ev, tracksH, assocOrigTrkRaw, assocOrigTrkToken);
780 }
781 
782 namespace {
783  bool cmp_for_detset(const unsigned one, const unsigned two) { return one < two; };
784 
785  void find_hits_in_dets(const MTDTrackingDetSetVector& hits,
786  const Trajectory& traj,
787  const DetLayer* layer,
788  const TrajectoryStateOnSurface& tsos,
789  const double pmag2,
790  const double pathlength0,
791  const double vtxTime,
792  const reco::BeamSpot& bs,
793  const float bsTimeSpread,
794  const Propagator* prop,
796  bool useVtxConstraint,
797  std::set<MTDHitMatchingInfo>& out) {
798  pair<bool, TrajectoryStateOnSurface> comp = layer->compatible(tsos, *prop, *estimator);
799  if (comp.first) {
800  const vector<DetLayer::DetWithState> compDets = layer->compatibleDets(tsos, *prop, *estimator);
801  if (!compDets.empty()) {
802  for (const auto& detWithState : compDets) {
803  auto range = hits.equal_range(detWithState.first->geographicalId(), cmp_for_detset);
804  if (range.first == range.second)
805  continue;
806 
807  auto pl = prop->propagateWithPath(tsos, detWithState.second.surface());
808  if (pl.second == 0.)
809  continue;
810 
811  const double tot_pl = pathlength0 + std::abs(pl.second);
812  const double t_vtx = useVtxConstraint ? vtxTime : 0.;
813 
814  constexpr double vtx_res = 0.008;
815  const double t_vtx_err = useVtxConstraint ? vtx_res : bsTimeSpread;
816 
817  constexpr double t_res_manual = 0.035;
818 
819  for (auto detitr = range.first; detitr != range.second; ++detitr) {
820  for (const auto& hit : *detitr) {
821  auto est = estimator->estimate(detWithState.second, hit);
822  if (!est.first)
823  continue;
824 
825  TrackTofPidInfo tof = computeTrackTofPidInfo(pmag2,
826  tot_pl,
827  hit.time(),
828  t_res_manual, //put hit error by hand for the moment
829  t_vtx,
830  t_vtx_err, //put vtx error by hand for the moment
831  false);
832  MTDHitMatchingInfo mi;
833  mi.hit = &hit;
834  mi.estChi2 = est.second;
835  mi.timeChi2 = tof.dtchi2_best; //use the chi2 for the best matching hypothesis
836 
837  out.insert(mi);
838  }
839  }
840  }
841  }
842  }
843  }
844 } // namespace
845 
846 template <class TrackCollection>
848  const TrajectoryStateOnSurface& tsos,
849  const Trajectory& traj,
850  const double pmag2,
851  const double pathlength0,
852  const MTDTrackingDetSetVector& hits,
853  const MTDDetLayerGeometry* geo,
854  const MagneticField* field,
855  const Propagator* prop,
856  const reco::BeamSpot& bs,
857  const double vtxTime,
858  const bool matchVertex,
859  MTDHitMatchingInfo& bestHit) const {
860  const vector<const DetLayer*>& layers = geo->allBTLLayers();
861 
863  bestHit = MTDHitMatchingInfo();
864  for (const DetLayer* ilay : layers)
865  fillMatchingHits(ilay, tsos, traj, pmag2, pathlength0, hits, prop, bs, vtxTime, matchVertex, output, bestHit);
866  return output;
867 }
868 
869 template <class TrackCollection>
871  const TrajectoryStateOnSurface& tsos,
872  const Trajectory& traj,
873  const double pmag2,
874  const double pathlength0,
875  const MTDTrackingDetSetVector& hits,
876  const MTDDetLayerGeometry* geo,
877  const MagneticField* field,
878  const Propagator* prop,
879  const reco::BeamSpot& bs,
880  const double vtxTime,
881  const bool matchVertex,
882  MTDHitMatchingInfo& bestHit) const {
883  const vector<const DetLayer*>& layers = geo->allETLLayers();
884 
886  bestHit = MTDHitMatchingInfo();
887  for (const DetLayer* ilay : layers) {
888  const BoundDisk& disk = static_cast<const ForwardDetLayer*>(ilay)->specificSurface();
889  const double diskZ = disk.position().z();
890 
891  if (tsos.globalPosition().z() * diskZ < 0)
892  continue; // only propagate to the disk that's on the same side
893 
894  fillMatchingHits(ilay, tsos, traj, pmag2, pathlength0, hits, prop, bs, vtxTime, matchVertex, output, bestHit);
895  }
896 
897  // the ETL hits order must be from the innermost to the outermost
898 
899  if (output.size() == 2) {
900  if (std::abs(output[0]->globalPosition().z()) > std::abs(output[1]->globalPosition().z())) {
901  std::reverse(output.begin(), output.end());
902  }
903  }
904  return output;
905 }
906 
907 template <class TrackCollection>
909  const TrajectoryStateOnSurface& tsos,
910  const Trajectory& traj,
911  const double pmag2,
912  const double pathlength0,
913  const MTDTrackingDetSetVector& hits,
914  const Propagator* prop,
915  const reco::BeamSpot& bs,
916  const double& vtxTime,
917  const bool matchVertex,
919  MTDHitMatchingInfo& bestHit) const {
920  std::set<MTDHitMatchingInfo> hitsInLayer;
921  bool hitMatched = false;
922 
923  using namespace std::placeholders;
924  auto find_hits = std::bind(find_hits_in_dets,
925  std::cref(hits),
926  std::cref(traj),
927  ilay,
928  std::cref(tsos),
929  pmag2,
930  pathlength0,
931  _1,
932  std::cref(bs),
933  bsTimeSpread_,
934  prop,
935  theEstimator.get(),
936  _2,
937  std::ref(hitsInLayer));
938 
939  if (useVertex_ && matchVertex)
940  find_hits(vtxTime, true);
941  else
942  find_hits(0, false);
943 
944  //just take the first hit because the hits are sorted on their matching quality
945  if (!hitsInLayer.empty()) {
946  //check hits to pass minimum quality matching requirements
947  auto const& firstHit = *hitsInLayer.begin();
948  if (firstHit.estChi2 < etlChi2Cut_ && firstHit.timeChi2 < etlTimeChi2Cut_) {
949  hitMatched = true;
950  output.push_back(hitbuilder_->build(firstHit.hit));
951  if (firstHit < bestHit)
952  bestHit = firstHit;
953  }
954  }
955 
956  if (useVertex_ && matchVertex && !hitMatched) {
957  //try a second search with beamspot hypothesis
958  hitsInLayer.clear();
959  find_hits(0, false);
960  if (!hitsInLayer.empty()) {
961  auto const& firstHit = *hitsInLayer.begin();
962  if (firstHit.timeChi2 < etlTimeChi2Cut_) {
963  if (firstHit.estChi2 < etlChi2Cut_) {
964  hitMatched = true;
965  output.push_back(hitbuilder_->build(firstHit.hit));
966  if (firstHit < bestHit)
967  bestHit = firstHit;
968  }
969  }
970  }
971  }
972 }
973 
974 //below is unfortunately ripped from other places but
975 //since track producer doesn't know about MTD we have to do this
976 template <class TrackCollection>
978  const Trajectory& traj,
979  const Trajectory& trajWithMtd,
980  const reco::BeamSpot& bs,
981  const MagneticField* field,
982  const Propagator* thePropagator,
983  bool hasMTD,
984  float& pathLengthOut,
985  float& tmtdOut,
986  float& sigmatmtdOut) const {
988  bool tsbcl_status = getTrajectoryStateClosestToBeamLine(traj, bs, thePropagator, tscbl);
989 
990  if (!tsbcl_status)
991  return reco::Track();
992 
994  math::XYZPoint pos(v.x(), v.y(), v.z());
996  math::XYZVector mom(p.x(), p.y(), p.z());
997 
998  int ndof = traj.ndof();
999 
1000  double t0 = 0.;
1001  double covt0t0 = -1.;
1002  pathLengthOut = -1.f; // if there is no MTD flag the pathlength with -1
1003  tmtdOut = 0.f;
1004  sigmatmtdOut = -1.f;
1005  double betaOut = 0.;
1006  double covbetabeta = -1.;
1007 
1008  auto routput = [&]() {
1009  return reco::Track(traj.chiSquared(),
1010  int(ndof),
1011  pos,
1012  mom,
1013  tscbl.trackStateAtPCA().charge(),
1015  orig.algo(),
1017  t0,
1018  betaOut,
1019  covt0t0,
1020  covbetabeta);
1021  };
1022 
1023  //compute path length for time backpropagation, using first MTD hit for the momentum
1024  if (hasMTD) {
1025  double pathlength;
1026  bool validpropagation = trackPathLength(trajWithMtd, bs, thePropagator, pathlength);
1027  double thit = 0.;
1028  double thiterror = -1.;
1029  bool validmtd = false;
1030 
1031  if (!validpropagation) {
1032  return routput();
1033  }
1034 
1035  size_t ihitcount(0), ietlcount(0);
1036  for (auto const& hit : trajWithMtd.measurements()) {
1037  if (hit.recHit()->geographicalId().det() == DetId::Forward &&
1038  ForwardSubdetector(hit.recHit()->geographicalId().subdetId()) == FastTime) {
1039  ihitcount++;
1040  if (MTDDetId(hit.recHit()->geographicalId()).mtdSubDetector() == MTDDetId::MTDType::ETL) {
1041  ietlcount++;
1042  }
1043  }
1044  }
1045 
1046  auto ihit1 = trajWithMtd.measurements().cbegin();
1047  if (ihitcount == 1) {
1048  const MTDTrackingRecHit* mtdhit = static_cast<const MTDTrackingRecHit*>((*ihit1).recHit()->hit());
1049  thit = mtdhit->time();
1050  thiterror = mtdhit->timeError();
1051  validmtd = true;
1052  } else if (ihitcount == 2 && ietlcount == 2) {
1053  const auto& propresult =
1054  thePropagator->propagateWithPath(ihit1->updatedState(), (ihit1 + 1)->updatedState().surface());
1055  double etlpathlength = std::abs(propresult.second);
1056  //
1057  // The information of the two ETL hits is combined and attributed to the innermost hit
1058  //
1059  if (etlpathlength == 0.) {
1060  validpropagation = false;
1061  } else {
1062  pathlength -= etlpathlength;
1063  const MTDTrackingRecHit* mtdhit1 = static_cast<const MTDTrackingRecHit*>((*ihit1).recHit()->hit());
1064  const MTDTrackingRecHit* mtdhit2 = static_cast<const MTDTrackingRecHit*>((*(ihit1 + 1)).recHit()->hit());
1065  TrackTofPidInfo tofInfo =
1066  computeTrackTofPidInfo(p.mag2(), etlpathlength, mtdhit1->time(), mtdhit1->timeError(), 0., 0., true);
1067  //
1068  // Protect against incompatible times
1069  //
1070  double err1 = tofInfo.dterror * tofInfo.dterror;
1071  double err2 = mtdhit2->timeError() * mtdhit2->timeError();
1072  if (cms_rounding::roundIfNear0(err1) == 0. || cms_rounding::roundIfNear0(err2) == 0.) {
1073  edm::LogError("TrackExtenderWithMTD")
1074  << "MTD tracking hits with zero time uncertainty: " << err1 << " " << err2;
1075  } else {
1076  if ((tofInfo.dt - mtdhit2->time()) * (tofInfo.dt - mtdhit2->time()) < (err1 + err2) * etlTimeChi2Cut_) {
1077  //
1078  // Subtract the ETL time of flight from the outermost measurement, and combine it in a weighted average with the innermost
1079  // the mass ambiguity related uncertainty on the time of flight is added as an additional uncertainty
1080  //
1081  err1 = 1. / err1;
1082  err2 = 1. / err2;
1083  thiterror = 1. / (err1 + err2);
1084  thit = (tofInfo.dt * err1 + mtdhit2->time() * err2) * thiterror;
1085  thiterror = std::sqrt(thiterror);
1086  LogDebug("TrackExtenderWithMTD") << "p trk = " << p.mag() << " ETL hits times/errors: " << mtdhit1->time()
1087  << " +/- " << mtdhit1->timeError() << " , " << mtdhit2->time() << " +/- "
1088  << mtdhit2->timeError() << " extrapolated time1: " << tofInfo.dt << " +/- "
1089  << tofInfo.dterror << " average = " << thit << " +/- " << thiterror;
1090  validmtd = true;
1091  }
1092  }
1093  }
1094  } else {
1095  edm::LogInfo("TrackExtenderWithMTD")
1096  << "MTD hits #" << ihitcount << "ETL hits #" << ietlcount << " anomalous pattern, skipping...";
1097  }
1098 
1099  if (validmtd && validpropagation) {
1100  //here add the PID uncertainty for later use in the 1st step of 4D vtx reconstruction
1101  TrackTofPidInfo tofInfo = computeTrackTofPidInfo(p.mag2(), pathlength, thit, thiterror, 0., 0., true);
1102  pathLengthOut = pathlength; // set path length if we've got a timing hit
1103  tmtdOut = thit;
1104  sigmatmtdOut = thiterror;
1105  t0 = tofInfo.dt;
1106  covt0t0 = tofInfo.dterror * tofInfo.dterror;
1107  betaOut = tofInfo.beta_pi;
1108  covbetabeta = tofInfo.betaerror * tofInfo.betaerror;
1109  }
1110  }
1111 
1112  return routput();
1113 }
1114 
1115 template <class TrackCollection>
1117  static const string metname = "TrackExtenderWithMTD";
1118 
1119  const Trajectory::RecHitContainer transRecHits = trajectory.recHits();
1120 
1121  // put the collection of TrackingRecHit in the event
1122 
1123  // sets the outermost and innermost TSOSs
1124  // ToDo: validation for track states with MTD
1125  TrajectoryStateOnSurface outerTSOS;
1126  TrajectoryStateOnSurface innerTSOS;
1127  unsigned int innerId = 0, outerId = 0;
1129  DetId outerDetId;
1130 
1131  if (trajectory.direction() == alongMomentum) {
1132  LogTrace(metname) << "alongMomentum";
1133  outerTSOS = trajectory.lastMeasurement().updatedState();
1134  innerTSOS = trajectory.firstMeasurement().updatedState();
1135  outerId = trajectory.lastMeasurement().recHit()->geographicalId().rawId();
1136  innerId = trajectory.firstMeasurement().recHit()->geographicalId().rawId();
1137  outerRecHit = trajectory.lastMeasurement().recHit();
1138  outerDetId = trajectory.lastMeasurement().recHit()->geographicalId();
1139  } else if (trajectory.direction() == oppositeToMomentum) {
1140  LogTrace(metname) << "oppositeToMomentum";
1141  outerTSOS = trajectory.firstMeasurement().updatedState();
1142  innerTSOS = trajectory.lastMeasurement().updatedState();
1143  outerId = trajectory.firstMeasurement().recHit()->geographicalId().rawId();
1144  innerId = trajectory.lastMeasurement().recHit()->geographicalId().rawId();
1145  outerRecHit = trajectory.firstMeasurement().recHit();
1146  outerDetId = trajectory.firstMeasurement().recHit()->geographicalId();
1147  } else
1148  LogError(metname) << "Wrong propagation direction!";
1149 
1150  const GeomDet* outerDet = gtg_->idToDet(outerDetId);
1151  GlobalPoint outerTSOSPos = outerTSOS.globalParameters().position();
1152  bool inside = outerDet->surface().bounds().inside(outerDet->toLocal(outerTSOSPos));
1153 
1154  GlobalPoint hitPos =
1155  (outerRecHit->isValid()) ? outerRecHit->globalPosition() : outerTSOS.globalParameters().position();
1156 
1157  if (!inside) {
1158  LogTrace(metname) << "The Global Muon outerMostMeasurementState is not compatible with the recHit detector!"
1159  << " Setting outerMost postition to recHit position if recHit isValid: "
1160  << outerRecHit->isValid();
1161  LogTrace(metname) << "From " << outerTSOSPos << " to " << hitPos;
1162  }
1163 
1164  //build the TrackExtra
1165  GlobalPoint v = (inside) ? outerTSOSPos : hitPos;
1166  GlobalVector p = outerTSOS.globalParameters().momentum();
1167  math::XYZPoint outpos(v.x(), v.y(), v.z());
1168  math::XYZVector outmom(p.x(), p.y(), p.z());
1169 
1170  v = innerTSOS.globalParameters().position();
1171  p = innerTSOS.globalParameters().momentum();
1172  math::XYZPoint inpos(v.x(), v.y(), v.z());
1173  math::XYZVector inmom(p.x(), p.y(), p.z());
1174 
1175  reco::TrackExtra trackExtra(outpos,
1176  outmom,
1177  true,
1178  inpos,
1179  inmom,
1180  true,
1181  outerTSOS.curvilinearError(),
1182  outerId,
1183  innerTSOS.curvilinearError(),
1184  innerId,
1185  trajectory.direction(),
1186  trajectory.seedRef());
1187 
1188  return trackExtra;
1189 }
1190 
1191 template <class TrackCollection>
1193  stringstream output;
1194 
1195  const BoundSurface* sur = nullptr;
1196  const BoundCylinder* bc = nullptr;
1197  const BoundDisk* bd = nullptr;
1198 
1199  sur = &(layer->surface());
1200  if ((bc = dynamic_cast<const BoundCylinder*>(sur))) {
1201  output << " Cylinder of radius: " << bc->radius() << endl;
1202  } else if ((bd = dynamic_cast<const BoundDisk*>(sur))) {
1203  output << " Disk at: " << bd->position().z() << endl;
1204  }
1205  return output.str();
1206 }
1207 
1208 //define this as a plug-in
1215 
1216 DEFINE_FWK_MODULE(TrackExtenderWithMTD);
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual std::pair< bool, TrajectoryStateOnSurface > compatible(const TrajectoryStateOnSurface &ts, const Propagator &, const MeasurementEstimator &) const =0
double z0() const
z coordinate
Definition: BeamSpot.h:65
float dt
Definition: AMPTWrapper.h:136
virtual const BoundSurface & surface() const =0
The surface of the GeometricSearchDet.
TrackExtenderWithMTDT< reco::TrackCollection > TrackExtenderWithMTD
T mag2() const
Definition: PV3DBase.h:63
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
edm::ESGetToken< TransientTrackingRecHitBuilder, TransientRecHitRecord > hitbuilderToken_
ConstRecHitPointer const & recHit() const
float time() const
Range equal_range(id_type i, CMP cmp, bool update=false) const
std::vector< unsigned char > Chi2sFive
reco::TrackExtra buildTrackExtra(const Trajectory &trajectory) const
constexpr double c_cm_ns
virtual bool inside(const Local3DPoint &) const =0
Determine if the point is inside the bounds.
edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > gtgToken_
const std::string metname
edm::EDGetTokenT< reco::BeamSpot > bsToken_
std::unique_ptr< MeasurementEstimator > theEstimator
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
TrackExtenderWithMTDT(const ParameterSet &pset)
const CurvilinearTrajectoryError & curvilinearError() const
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > ttopoToken_
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
edm::EDPutToken pathLengthOrigTrkToken
T y() const
Definition: PV3DBase.h:60
string dumpLayer(const DetLayer *layer) const
edm::View< TrackType > InputCollection
edm::EDGetTokenT< GlobalPoint > genVtxPositionToken_
GlobalPoint globalPosition() const
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
ConstRecHitContainer recHits() const
Definition: Trajectory.h:186
const Bounds & bounds() const
Definition: Surface.h:87
auto const & tracks
cannot be loose
bool ev
const std::vector< const DetLayer * > & allBTLLayers() const
return the BTL DetLayers (barrel), inside-out
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
TrackCharge charge() const
reco::Track buildTrack(const reco::Track &, const Trajectory &, const Trajectory &, const reco::BeamSpot &, const MagneticField *field, const Propagator *prop, bool hasMTD, float &pathLength, float &tmtdOut, float &sigmatmtdOut) const
Log< level::Error, false > LogError
TrajectoryMeasurement const & closestMeasurement(GlobalPoint) const
Definition: Trajectory.cc:256
ForwardSubdetector
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
const Point & position() const
position
Definition: Vertex.h:127
const CurvilinearTrajectoryError & curvilinearError() const
float float float z
virtual std::vector< DetWithState > compatibleDets(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
edm::EDPutToken assocOrigTrkToken
edm::EDGetTokenT< MTDTrackingDetSetVector > hitsToken_
#define LogTrace(id)
Detector identifier base class for the MIP Timing Layer.
Definition: MTDDetId.h:21
constexpr std::array< uint8_t, layerIndexSize > layer
tuple m2
Definition: callgraph.py:57
const uint16_t range(const Frame &aFrame)
tuple result
Definition: mps_fire.py:311
PropagationDirection const & direction() const
Definition: Trajectory.cc:133
TrackAlgorithm algo() const
Definition: TrackBase.h:547
list ordering
Definition: config.py:7
DataContainer const & measurements() const
Definition: Trajectory.h:178
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:563
std::vector< LocalTrajectoryParameters > TrajParams
int iEvent
Definition: GenABIO.cc:224
const SurfaceType & surface() const
edm::EDPutToken tmtdOrigTrkToken
TransientTrackingRecHit::ConstRecHitContainer tryETLLayers(const TrajectoryStateOnSurface &, const Trajectory &traj, const double, const double, const MTDTrackingDetSetVector &, const MTDDetLayerGeometry *, const MagneticField *field, const Propagator *prop, const reco::BeamSpot &bs, const double vtxTime, const bool matchVertex, MTDHitMatchingInfo &bestHit) const
edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > builderToken_
const std::vector< const DetLayer * > & allETLLayers() const
return the ETL DetLayers (endcap), -Z to +Z
const std::string mtdRecHitBuilder_
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:590
T sqrt(T t)
Definition: SSEVec.h:19
edm::EDPutToken betaOrigTrkToken
TrajectoryMeasurement const & lastMeasurement() const
Definition: Trajectory.h:150
FreeTrajectoryState const * freeState(bool withErrors=true) const
T z() const
Definition: PV3DBase.h:61
def move
Definition: eostools.py:511
RefitDirection::GeometricalDirection checkRecHitsOrdering(TransientTrackingRecHit::ConstRecHitContainer const &recHits) const
tuple handle
Definition: patZpeak.py:23
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:346
const std::string propagator_
edm::ESGetToken< Propagator, TrackingComponentsRecord > propToken_
virtual HitReturnType estimate(const TrajectoryStateOnSurface &ts, const TrackingRecHit &hit) const =0
edm::EDPutToken etlMatchChi2Token
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::RefToBase< TrajectorySeed > seedRef(void) const
Definition: Trajectory.h:303
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
edm::ESHandle< TransientTrackingRecHitBuilder > hitbuilder_
GlobalVector momentum() const
RefProd< PROD > getRefBeforePut()
Definition: Event.h:158
std::vector< ConstRecHitPointer > ConstRecHitContainer
constexpr double c_inv
const std::string transientTrackBuilder_
Log< level::Info, false > LogInfo
Definition: DetId.h:17
GlobalPoint position() const
int ndof(bool bon=true) const
Definition: Trajectory.cc:97
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:10
TrajectoryMeasurement const & firstMeasurement() const
Definition: Trajectory.h:166
ConstRecHitContainer RecHitContainer
Definition: Trajectory.h:42
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const GlobalTrajectoryParameters & globalParameters() const
ParameterSet const & getParameterSet(std::string const &) const
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::EDPutToken sigmatmtdOrigTrkToken
A 2D TrackerRecHit with time and time error information.
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void add(std::string const &label, ParameterSetDescription const &psetDescription)
TrackingRecHit::ConstRecHitPointer ConstRecHitPointer
edm::ESHandle< TransientTrackBuilder > builder_
TransientTrackingRecHit::ConstRecHitContainer tryBTLLayers(const TrajectoryStateOnSurface &, const Trajectory &traj, const double, const double, const MTDTrackingDetSetVector &, const MTDDetLayerGeometry *, const MagneticField *field, const Propagator *prop, const reco::BeamSpot &bs, const double vtxTime, const bool matchVertex, MTDHitMatchingInfo &bestHit) const
edm::EDPutToken btlMatchTimeChi2Token
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
edm::EDPutToken etlMatchTimeChi2Token
std::unique_ptr< TrackTransformer > theTransformer
float chiSquared() const
Definition: Trajectory.h:241
edm::EDPutToken sigmat0OrigTrkToken
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magfldToken_
bool operator<(DTCELinkId const &lhs, DTCELinkId const &rhs)
Definition: DTCELinkId.h:70
static int position[264][3]
Definition: ReadPGInfo.cc:289
ESTransientHandle< T > getTransientHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:168
constexpr NumType convertMmToCm(NumType millimeters)
Definition: GeantUnits.h:63
double y0() const
y coordinate
Definition: BeamSpot.h:63
TrajectoryStateOnSurface const & updatedState() const
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:157
tuple last
Definition: dqmdumpme.py:56
edm::ESHandle< GlobalTrackingGeometry > gtg_
float timeError() const
edm::EDGetTokenT< float > genVtxTimeToken_
T x() const
Definition: PV3DBase.h:59
edm::EDGetTokenT< VertexCollection > vtxToken_
edm::EDPutToken btlMatchChi2Token
void produce(edm::Event &ev, const edm::EventSetup &es) final
TrackCollection::value_type TrackType
edm::EDGetTokenT< InputCollection > tracksToken_
double t() const
t coordinate
Definition: Vertex.h:135
void fillValueMap(edm::Event &iEvent, const H &handle, const std::vector< T > &vec, const edm::EDPutToken &token) const
void fillMatchingHits(const DetLayer *, const TrajectoryStateOnSurface &, const Trajectory &, const double, const double, const MTDTrackingDetSetVector &, const Propagator *, const reco::BeamSpot &, const double &, const bool, TransientTrackingRecHit::ConstRecHitContainer &, MTDHitMatchingInfo &) const
#define LogDebug(id)
edm::ESGetToken< MTDDetLayerGeometry, MTDRecoGeometryRecord > dlgeoToken_
#define m_pi
Definition: RPCConst.cc:8
double x0() const
x coordinate
Definition: BeamSpot.h:61