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