CMS 3D CMS Logo

TrackTransformerForCosmicMuons.cc
Go to the documentation of this file.
2 
6 
9 
13 
20 
24 
29 
30 using namespace std;
31 using namespace edm;
32 
35  theTrackerRecHitBuilderName = parameterSet.getParameter<string>("TrackerRecHitBuilder");
36  theMuonRecHitBuilderName = parameterSet.getParameter<string>("MuonRecHitBuilder");
37 
38  theRPCInTheFit = parameterSet.getParameter<bool>("RefitRPCHits");
39 
40  theCacheId_TC = theCacheId_GTG = theCacheId_MG = theCacheId_TRH = 0;
41 }
42 
45 
47  const std::string metname = "Reco|TrackingTools|TrackTransformer";
48 
49  setup.get<TrajectoryFitter::Record>().get("KFFitterForRefitInsideOut", theFitterIO);
50  setup.get<TrajectoryFitter::Record>().get("KFSmootherForRefitInsideOut", theSmootherIO);
51  setup.get<TrajectoryFitter::Record>().get("KFFitterForRefitOutsideIn", theFitterOI);
52  setup.get<TrajectoryFitter::Record>().get("KFSmootherForRefitOutsideIn", theSmootherOI);
53 
54  unsigned long long newCacheId_TC = setup.get<TrackingComponentsRecord>().cacheIdentifier();
55 
56  if (newCacheId_TC != theCacheId_TC) {
57  LogTrace(metname) << "Tracking Component changed!";
58  theCacheId_TC = newCacheId_TC;
59  setup.get<TrackingComponentsRecord>().get("SmartPropagatorRK", thePropagatorIO);
60  setup.get<TrackingComponentsRecord>().get("SmartPropagatorRKOpposite", thePropagatorOI);
61  }
62 
63  // Global Tracking Geometry
64  unsigned long long newCacheId_GTG = setup.get<GlobalTrackingGeometryRecord>().cacheIdentifier();
65  if (newCacheId_GTG != theCacheId_GTG) {
66  LogTrace(metname) << "GlobalTrackingGeometry changed!";
67  theCacheId_GTG = newCacheId_GTG;
68  setup.get<GlobalTrackingGeometryRecord>().get(theTrackingGeometry);
69  }
70 
71  // Magfield Field
72  unsigned long long newCacheId_MG = setup.get<IdealMagneticFieldRecord>().cacheIdentifier();
73  if (newCacheId_MG != theCacheId_MG) {
74  LogTrace(metname) << "Magnetic Field changed!";
75  theCacheId_MG = newCacheId_MG;
76  setup.get<IdealMagneticFieldRecord>().get(theMGField);
77  }
78 
79  // Transient Rechit Builders
80  unsigned long long newCacheId_TRH = setup.get<TransientRecHitRecord>().cacheIdentifier();
81  if (newCacheId_TRH != theCacheId_TRH) {
82  theCacheId_TRH = newCacheId_TRH;
83  LogTrace(metname) << "TransientRecHitRecord changed!";
84  setup.get<TransientRecHitRecord>().get(theTrackerRecHitBuilderName, theTrackerRecHitBuilder);
85  setup.get<TransientRecHitRecord>().get(theMuonRecHitBuilderName, theMuonRecHitBuilder);
86  }
87 }
88 
90  const reco::TransientTrack& track) const {
93 
94  bool printout = false;
95 
96  bool quad1 = false;
97  bool quad2 = false;
98  bool quad3 = false;
99  bool quad4 = false;
100 
101  for (trackingRecHit_iterator hit = track.recHitsBegin(); hit != track.recHitsEnd(); ++hit)
102  if ((*hit)->isValid())
103  if ((*hit)->geographicalId().det() == DetId::Muon) {
104  if ((*hit)->geographicalId().subdetId() == 3 && !theRPCInTheFit) {
105  LogTrace("Reco|TrackingTools|TrackTransformer") << "RPC Rec Hit discarged";
106  continue;
107  }
108  staHits.push_back(theMuonRecHitBuilder->build(&**hit));
109  DetId hitId = staHits.back()->geographicalId();
110  GlobalPoint glbpoint = trackingGeometry()->idToDet(hitId)->position();
111  float gpy = glbpoint.y();
112  float gpz = glbpoint.z();
113  // if (gpy != 0 && gpz !=0) slopeSum += gpy / gpz;
114 
115  if (gpy > 0 && gpz > 0)
116  quad1 = true;
117  else if (gpy > 0 && gpz < 0)
118  quad2 = true;
119  else if (gpy < 0 && gpz < 0)
120  quad3 = true;
121  else if (gpy < 0 && gpz > 0)
122  quad4 = true;
123  else
124  return tkHits;
125  }
126 
127  if (staHits.empty())
128  return staHits;
129 
130  if (quad1 && quad2)
131  return tkHits;
132  if (quad1 && quad3)
133  return tkHits;
134  if (quad1 && quad4)
135  return tkHits;
136  if (quad2 && quad3)
137  return tkHits;
138  if (quad2 && quad4)
139  return tkHits;
140  if (quad3 && quad4)
141  return tkHits;
142 
143  bool doReverse = staHits.front()->globalPosition().y() > 0 ? true : false;
144 
145  if (quad1)
146  doReverse = SlopeSum(staHits);
147  if (quad2)
148  doReverse = SlopeSum(staHits);
149  if (quad3)
150  doReverse = SlopeSum(staHits);
151  if (quad4)
152  doReverse = SlopeSum(staHits);
153  if (doReverse) {
154  reverse(staHits.begin(), staHits.end());
155  }
156 
157  copy(staHits.begin(), staHits.end(), back_inserter(tkHits));
158 
160  // if ( quad2 && slopeSum > 0) printout = true;
163  // swap = printout;
164 
165  printout = quad1;
166 
167  if (printout)
168  for (TransientTrackingRecHit::ConstRecHitContainer::const_iterator hit = tkHits.begin(); hit != tkHits.end();
169  ++hit) {
170  DetId hitId = (*hit)->geographicalId();
171  GlobalPoint glbpoint = trackingGeometry()->idToDet(hitId)->position();
172 
173  if (hitId.det() == DetId::Muon) {
174  if (hitId.subdetId() == MuonSubdetId::DT)
175  LogTrace("TrackFitters") << glbpoint << " I am DT " << DTWireId(hitId);
176  // std::cout<< glbpoint << " I am DT " << DTWireId(hitId)<<std::endl;
177  else if (hitId.subdetId() == MuonSubdetId::CSC)
178  LogTrace("TrackFitters") << glbpoint << " I am CSC " << CSCDetId(hitId);
179  // std::cout<< glbpoint << " I am CSC " << CSCDetId(hitId)<<std::endl;
180  else if (hitId.subdetId() == MuonSubdetId::RPC)
181  LogTrace("TrackFitters") << glbpoint << " I am RPC " << RPCDetId(hitId);
182  else
183  LogTrace("TrackFitters") << " UNKNOWN MUON HIT TYPE ";
184  }
185  }
186  return tkHits;
187 }
188 
191  if (quad == 1) {
192  if (sumy < 0)
193  return theFitterOI;
194  else
195  return theFitterIO;
196  }
197  if (quad == 2) {
198  if (sumy < 0)
199  return theFitterOI;
200  else
201  return theFitterIO;
202  }
203  if (quad == 3) {
204  if (sumy > 0)
205  return theFitterOI;
206  else
207  return theFitterIO;
208  }
209  if (quad == 4) {
210  if (sumy > 0)
211  return theFitterOI;
212  else
213  return theFitterIO;
214  }
215 
216  if (up)
217  return theFitterOI;
218  else
219  return theFitterIO;
220 }
221 
224  if (quad == 1) {
225  if (sumy < 0)
226  return theSmootherOI;
227  else
228  return theSmootherIO;
229  }
230  if (quad == 2) {
231  if (sumy < 0)
232  return theSmootherOI;
233  else
234  return theSmootherIO;
235  }
236  if (quad == 3) {
237  if (sumy > 0)
238  return theSmootherOI;
239  else
240  return theSmootherIO;
241  }
242  if (quad == 4) {
243  if (sumy > 0)
244  return theSmootherOI;
245  else
246  return theSmootherIO;
247  }
248  if (up)
249  return theSmootherOI;
250  else
251  return theSmootherIO;
252 }
253 
255  if (quad == 1) {
256  if (sumy > 0)
257  return thePropagatorOI;
258  else
259  return thePropagatorIO;
260  }
261  if (quad == 2) {
262  if (sumy > 0)
263  return thePropagatorOI;
264  else
265  return thePropagatorIO;
266  }
267  if (quad == 3) {
268  if (sumy < 0)
269  return thePropagatorOI;
270  else
271  return thePropagatorIO;
272  }
273  if (quad == 4) {
274  if (sumy < 0)
275  return thePropagatorOI;
276  else
277  return thePropagatorIO;
278  }
279  if (up)
280  return thePropagatorIO;
281  else
282  return thePropagatorOI;
283 }
284 
286 vector<Trajectory> TrackTransformerForCosmicMuons::transform(const reco::Track& tr) const {
287  const std::string metname = "Reco|TrackingTools|TrackTransformer";
288 
289  reco::TransientTrack track(tr, magneticField(), trackingGeometry());
290 
291  // Build the transient Rechits
292  TransientTrackingRecHit::ConstRecHitContainer recHitsForReFit; // = getTransientRecHits(track);
294 
295  bool quad1 = false;
296  bool quad2 = false;
297  bool quad3 = false;
298  bool quad4 = false;
299  int quadrant = 0;
300 
301  for (trackingRecHit_iterator hit = track.recHitsBegin(); hit != track.recHitsEnd(); ++hit)
302  if ((*hit)->isValid())
303  if ((*hit)->geographicalId().det() == DetId::Muon) {
304  if ((*hit)->geographicalId().subdetId() == 3 && !theRPCInTheFit) {
305  LogTrace("Reco|TrackingTools|TrackTransformer") << "RPC Rec Hit discarged";
306  continue;
307  }
308  staHits.push_back(theMuonRecHitBuilder->build(&**hit));
309  DetId hitId = staHits.back()->geographicalId();
310  GlobalPoint glbpoint = trackingGeometry()->idToDet(hitId)->position();
311  float gpy = glbpoint.y();
312  float gpz = glbpoint.z();
313  // if (gpy != 0 && gpz !=0) slopeSum += gpy / gpz;
314 
315  if (gpy > 0 && gpz > 0) {
316  quad1 = true;
317  quadrant = 1;
318  } else if (gpy > 0 && gpz < 0) {
319  quad2 = true;
320  quadrant = 2;
321  } else if (gpy < 0 && gpz < 0) {
322  quad3 = true;
323  quadrant = 3;
324  } else if (gpy < 0 && gpz > 0) {
325  quad4 = true;
326  quadrant = 4;
327  } else
328  return vector<Trajectory>();
329  }
330 
331  if (staHits.empty())
332  return vector<Trajectory>();
333 
334  if (quad1 && quad2)
335  return vector<Trajectory>();
336  if (quad1 && quad3)
337  return vector<Trajectory>();
338  if (quad1 && quad4)
339  return vector<Trajectory>();
340  if (quad2 && quad3)
341  return vector<Trajectory>();
342  if (quad2 && quad4)
343  return vector<Trajectory>();
344  if (quad3 && quad4)
345  return vector<Trajectory>();
346 
347  bool doReverse = false;
348 
349  if (quad1)
350  doReverse = SlopeSum(staHits);
351  if (quad2)
352  doReverse = SlopeSum(staHits);
353  if (quad3)
354  doReverse = SlopeSum(staHits);
355  if (quad4)
356  doReverse = SlopeSum(staHits);
357 
358  if (doReverse) {
359  reverse(staHits.begin(), staHits.end());
360  }
361 
362  copy(staHits.begin(), staHits.end(), back_inserter(recHitsForReFit));
363 
369 
370  if (recHitsForReFit.size() < 2)
371  return vector<Trajectory>();
372 
373  bool up = recHitsForReFit.back()->globalPosition().y() > 0 ? true : false;
374  LogTrace(metname) << "Up ? " << up;
375 
376  float sumdy = SumDy(recHitsForReFit);
377 
378  PropagationDirection propagationDirection = doReverse ? oppositeToMomentum : alongMomentum;
379  TrajectoryStateOnSurface firstTSOS =
380  doReverse ? track.outermostMeasurementState() : track.innermostMeasurementState();
381  unsigned int innerId = doReverse ? track.track().outerDetId() : track.track().innerDetId();
382 
383  LogTrace(metname) << "Prop Dir: " << propagationDirection << " FirstId " << innerId << " firstTSOS " << firstTSOS;
384 
386 
387  if (recHitsForReFit.front()->geographicalId() != DetId(innerId)) {
388  LogTrace(metname) << "Propagation occurring" << endl;
389  firstTSOS = propagator(up, quadrant, sumdy)->propagate(firstTSOS, recHitsForReFit.front()->det()->surface());
390  LogTrace(metname) << "Final destination: " << recHitsForReFit.front()->det()->surface().position() << endl;
391  if (!firstTSOS.isValid()) {
392  std::cout << "Propagation error! Dumping RecHits..." << std::endl;
393 
394  for (TransientTrackingRecHit::ConstRecHitContainer::const_iterator hit = recHitsForReFit.begin();
395  hit != recHitsForReFit.end();
396  ++hit) {
397  DetId hitId = (*hit)->geographicalId();
398  GlobalPoint glbpoint = trackingGeometry()->idToDet(hitId)->position();
399  if (hitId.subdetId() == MuonSubdetId::DT)
400  std::cout << glbpoint << " I am DT " << DTWireId(hitId) << std::endl;
401  else if (hitId.subdetId() == MuonSubdetId::CSC)
402  std::cout << glbpoint << " I am CSC " << CSCDetId(hitId) << std::endl;
403  }
404 
405  LogTrace(metname) << "Propagation error!" << endl;
406  return vector<Trajectory>();
407  }
408  }
409 
410  vector<Trajectory>&& trajectories = fitter(up, quadrant, sumdy)->fit(seed, recHitsForReFit, firstTSOS);
411 
412  if (trajectories.empty()) {
413  LogTrace(metname) << "No Track refitted!" << endl;
414  return vector<Trajectory>();
415  }
416 
417  Trajectory const& trajectoryBW = trajectories.front();
418 
419  vector<Trajectory> trajectoriesSM = smoother(up, quadrant, sumdy)->trajectories(trajectoryBW);
420 
421  if (trajectoriesSM.empty()) {
422  LogTrace(metname) << "No Track smoothed!" << endl;
423  }
424 
425  return trajectoriesSM;
426 }
427 
430  bool retval = false;
431  float y1 = 0;
432  float z1 = 0;
433 
434  bool first = true;
435 
436  std::vector<float> py;
437  std::vector<float> pz;
438  // int quadrant = -1;
439 
440  float sumdy = 0;
441  float sumdz = 0;
442 
443  for (TransientTrackingRecHit::ConstRecHitContainer::const_iterator hit = tkHits.begin(); hit != tkHits.end(); ++hit) {
444  DetId hitId = (*hit)->geographicalId();
445  GlobalPoint glbpoint = trackingGeometry()->idToDet(hitId)->position();
446  if (hitId.det() != DetId::Muon || hitId.subdetId() == 3)
447  continue;
448 
449  float y2 = glbpoint.y();
450  float z2 = glbpoint.z();
451  float dslope = 0;
452  float dy = y2 - y1;
453  float dz = z2 - z1;
454 
455  // if (y2 > 0 && z2 > 0) quadrant = 1;
456  // else if (y2 > 0 && z2 < 0) quadrant = 2;
457  // else if (y2 < 0 && z2 < 0) quadrant = 3;
458  // else if (y2 < 0 && z2 > 0) quadrant = 4;
459 
460  if (fabs(dz) > 1e-3)
461  dslope = dy / dz;
462  if (!first) {
463  retval += dslope;
464  sumdy += dy;
465  sumdz += dz;
466  }
467  first = false;
468  py.push_back(y1);
469  pz.push_back(z1);
470  y1 = y2;
471  z1 = z2;
472  }
473  // std::cout<<"quadrant "<<quadrant;
474  // std::cout<<"\tsum dy = "<<sumdy;
475  // std::cout<<"\tsum dz = "<<sumdz;
476  // std::cout<<std::endl;
477 
478  if (sumdz < 0)
479  retval = true;
480 
481  return retval;
482 }
483 
486  bool retval = false;
487  float y1 = 0;
488  float z1 = 0;
489 
490  bool first = true;
491 
492  std::vector<float> py;
493  std::vector<float> pz;
494  // int quadrant = -1;
495 
496  float sumdy = 0;
497  float sumdz = 0;
498 
499  for (TransientTrackingRecHit::ConstRecHitContainer::const_iterator hit = tkHits.begin(); hit != tkHits.end(); ++hit) {
500  DetId hitId = (*hit)->geographicalId();
501  GlobalPoint glbpoint = trackingGeometry()->idToDet(hitId)->position();
502  if (hitId.det() != DetId::Muon || hitId.subdetId() == 3)
503  continue;
504 
505  float y2 = glbpoint.y();
506  float z2 = glbpoint.z();
507  float dslope = 0;
508  float dy = y2 - y1;
509  float dz = z2 - z1;
510 
511  // if (y2 > 0 && z2 > 0) quadrant = 1;
512  // else if (y2 > 0 && z2 < 0) quadrant = 2;
513  // else if (y2 < 0 && z2 < 0) quadrant = 3;
514  // else if (y2 < 0 && z2 > 0) quadrant = 4;
515 
516  if (fabs(dz) > 1e-3)
517  dslope = dy / dz;
518  if (!first) {
519  retval += dslope;
520  sumdy += dy;
521  sumdz += dz;
522  }
523  first = false;
524  py.push_back(y1);
525  pz.push_back(z1);
526  y1 = y2;
527  z1 = z2;
528  }
529  // std::cout<<"quadrant "<<quadrant;
530  // std::cout<<"\tsum dy = "<<sumdy;
531  // std::cout<<"\tsum dz = "<<sumdz;
532  // std::cout<<std::endl;
533 
534  return sumdy;
535 }
Definition: BitonicSort.h:7
T getParameter(std::string const &) const
std::vector< Trajectory > transform(const reco::Track &) const override
Convert a reco::Track into Trajectory.
edm::ESHandle< Propagator > propagator(bool, int, float) const
void setServices(const edm::EventSetup &) override
set the services needed by the TrackTransformer
const std::string metname
int quadrant(const DetId &detid, const TrackerTopology *tTopo_, bool phase_)
T y() const
Definition: PV3DBase.h:60
TrackTransformerForCosmicMuons(const edm::ParameterSet &)
Constructor.
PropagationDirection
edm::ESHandle< TrajectoryFitter > fitter(bool, int, float) const
the refitter used to refit the reco::Track
TrajectoryStateOnSurface innermostMeasurementState() const
float SumDy(const TransientTrackingRecHit::ConstRecHitContainer &) const
decide if the track should be reversed
~TrackTransformerForCosmicMuons() override
Destructor.
T z() const
Definition: PV3DBase.h:61
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
unsigned int outerDetId() const
DetId of the detector on which surface the outermost state is located.
Definition: Track.h:79
bool SlopeSum(const TransientTrackingRecHit::ConstRecHitContainer &) const
calculate the sum of slopes for the track
TrajectoryStateOnSurface outermostMeasurementState() const
#define LogTrace(id)
trackingRecHit_iterator recHitsEnd() const
last iterator to RecHits
std::vector< ConstRecHitPointer > ConstRecHitContainer
Definition: DetId.h:17
const Track & track() const
static constexpr int RPC
Definition: MuonSubdetId.h:13
edm::ESHandle< TrajectorySmoother > smoother(bool, int, float) const
the smoother used to smooth the trajectory which came from the refitting step
HLT enums.
T get() const
Definition: EventSetup.h:73
static constexpr int DT
Definition: MuonSubdetId.h:11
static constexpr int CSC
Definition: MuonSubdetId.h:12
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:82
trackingRecHit_iterator recHitsBegin() const
first iterator to RecHits
ParameterSet const & parameterSet(Provenance const &provenance)
Definition: Provenance.cc:11
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:46
TransientTrackingRecHit::ConstRecHitContainer getTransientRecHits(const reco::TransientTrack &track) const