CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GlobalMuonRefitter.cc
Go to the documentation of this file.
1 
16 
17 //---------------
18 // C++ Headers --
19 //---------------
20 
21 #include <iostream>
22 #include <iomanip>
23 #include <algorithm>
24 
25 //-------------------------------
26 // Collaborating Class Headers --
27 //-------------------------------
28 
30 
40 
41 
56 
59 
67 
68 using namespace std;
69 using namespace edm;
70 
71 //----------------
72 // Constructors --
73 //----------------
74 
76  const MuonServiceProxy* service) :
77  theCosmicFlag(par.getParameter<bool>("PropDirForCosmics")),
78  theDTRecHitLabel(par.getParameter<InputTag>("DTRecSegmentLabel")),
79  theCSCRecHitLabel(par.getParameter<InputTag>("CSCRecSegmentLabel")),
80  theService(service) {
81 
82  theCategory = par.getUntrackedParameter<string>("Category", "Muon|RecoMuon|GlobalMuon|GlobalMuonRefitter");
83 
84  theHitThreshold = par.getParameter<int>("HitThreshold");
85  theDTChi2Cut = par.getParameter<double>("Chi2CutDT");
86  theCSCChi2Cut = par.getParameter<double>("Chi2CutCSC");
87  theRPCChi2Cut = par.getParameter<double>("Chi2CutRPC");
88 
89  // Refit direction
90  string refitDirectionName = par.getParameter<string>("RefitDirection");
91 
92  if (refitDirectionName == "insideOut" ) theRefitDirection = insideOut;
93  else if (refitDirectionName == "outsideIn" ) theRefitDirection = outsideIn;
94  else
95  throw cms::Exception("TrackTransformer constructor")
96  <<"Wrong refit direction chosen in TrackTransformer ParameterSet"
97  << "\n"
98  << "Possible choices are:"
99  << "\n"
100  << "RefitDirection = insideOut or RefitDirection = outsideIn";
101 
102  theFitterName = par.getParameter<string>("Fitter");
103  thePropagatorName = par.getParameter<string>("Propagator");
104 
105  theSkipStation = par.getParameter<int>("SkipStation");
106  theTrackerSkipSystem = par.getParameter<int>("TrackerSkipSystem");
107  theTrackerSkipSection = par.getParameter<int>("TrackerSkipSection");//layer, wheel, or disk depending on the system
108 
109  theTrackerRecHitBuilderName = par.getParameter<string>("TrackerRecHitBuilder");
110  theMuonRecHitBuilderName = par.getParameter<string>("MuonRecHitBuilder");
111 
112  theRPCInTheFit = par.getParameter<bool>("RefitRPCHits");
113 
114  theCacheId_TRH = 0;
115 
116 }
117 
118 //--------------
119 // Destructor --
120 //--------------
121 
123 }
124 
125 
126 //
127 // set Event
128 //
130 
131  theEvent = &event;
132  event.getByLabel(theDTRecHitLabel, theDTRecHits);
133  event.getByLabel(theCSCRecHitLabel, theCSCRecHits);
134 }
135 
136 
138 
140 
141  // Transient Rechit Builders
142  unsigned long long newCacheId_TRH = setup.get<TransientRecHitRecord>().cacheIdentifier();
143  if ( newCacheId_TRH != theCacheId_TRH ) {
144  LogDebug(theCategory) << "TransientRecHitRecord changed!";
147  }
148 }
149 
150 
151 //
152 // build a combined tracker-muon trajectory
153 //
154 vector<Trajectory> GlobalMuonRefitter::refit(const reco::Track& globalTrack,
155  const int theMuonHitsOption) const {
156  LogTrace(theCategory) << " *** GlobalMuonRefitter *** option " << theMuonHitsOption << endl;
157 
158  ConstRecHitContainer allRecHitsTemp; // all muon rechits temp
159 
160  reco::TransientTrack track(globalTrack,&*(theService->magneticField()),theService->trackingGeometry());
161 
162  for (trackingRecHit_iterator hit = track.recHitsBegin(); hit != track.recHitsEnd(); ++hit)
163  if ((*hit)->isValid()) {
164  if ((*hit)->geographicalId().det() == DetId::Tracker)
165  allRecHitsTemp.push_back(theTrackerRecHitBuilder->build(&**hit));
166  else if ((*hit)->geographicalId().det() == DetId::Muon) {
167  if ((*hit)->geographicalId().subdetId() == 3 && !theRPCInTheFit) {
168  LogTrace(theCategory) << "RPC Rec Hit discarged";
169  continue;
170  }
171  allRecHitsTemp.push_back(theMuonRecHitBuilder->build(&**hit));
172  }
173  }
174  vector<Trajectory> refitted = refit(globalTrack,track,allRecHitsTemp,theMuonHitsOption);
175  return refitted;
176 }
177 
178 //
179 // build a combined tracker-muon trajectory
180 //
181 vector<Trajectory> GlobalMuonRefitter::refit(const reco::Track& globalTrack,
182  const reco::TransientTrack track,
184  const int theMuonHitsOption) const {
185 
186  // MuonHitsOption: 0 - tracker only
187  // 1 - include all muon hits
188  // 2 - include only first muon hit(s)
189  // 3 - include only selected muon hits
190  // 4 - redo pattern recognition with dynamic truncation
191 
192  vector<int> stationHits(4,0);
193 
194  ConstRecHitContainer allRecHits; // all muon rechits
195  ConstRecHitContainer fmsRecHits; // only first muon rechits
196  ConstRecHitContainer selectedRecHits; // selected muon rechits
197  ConstRecHitContainer DYTRecHits; // rec hits from dynamic truncation algorithm
198 
199  LogTrace(theCategory) << " *** GlobalMuonRefitter *** option " << theMuonHitsOption << endl;
200 
201  LogTrace(theCategory) << " Track momentum before refit: " << globalTrack.pt() << endl;
202 
203  LogTrace(theCategory) << " Hits size before : " << allRecHitsTemp.size() << endl;
204  allRecHits = getRidOfSelectStationHits(allRecHitsTemp);
205  // printHits(allRecHits);
206  LogTrace(theCategory) << " Hits size: " << allRecHits.size() << endl;
207 
208  vector <Trajectory> outputTraj;
209 
210  if ((theMuonHitsOption == 1) || (theMuonHitsOption == 3) || (theMuonHitsOption == 4) ) {
211  // refit the full track with all muon hits
212  vector <Trajectory> globalTraj = transform(globalTrack, track, allRecHits);
213 
214  if (!globalTraj.size()) {
215  LogTrace(theCategory) << "No trajectory from the TrackTransformer!" << endl;
216  return vector<Trajectory>();
217  }
218 
219  LogTrace(theCategory) << " Initial trajectory state: "
220  << globalTraj.front().lastMeasurement().updatedState().freeState()->parameters() << endl;
221 
222  if (theMuonHitsOption == 1 )
223  outputTraj.push_back(globalTraj.front());
224 
225  if (theMuonHitsOption == 3 ) {
226  checkMuonHits(globalTrack, allRecHits, stationHits);
227  selectedRecHits = selectMuonHits(globalTraj.front(),stationHits);
228  LogTrace(theCategory) << " Selected hits size: " << selectedRecHits.size() << endl;
229  outputTraj = transform(globalTrack, track, selectedRecHits);
230  }
231 
232  if (theMuonHitsOption == 4 ) {
233  // here we use the single thr per subdetector (better performance can be obtained using thr as function of eta)
234 
236  dytRefit.setThr(30,15);
237  //dytRefit.setThr(20,20,20,20,20,15,15,15,15,15,15,15,15);
238  DYTRecHits = dytRefit.filter(globalTraj.front());
239  if ((DYTRecHits.size() > 1) && (DYTRecHits.front()->globalPosition().mag() > DYTRecHits.back()->globalPosition().mag()))
240  stable_sort(DYTRecHits.begin(),DYTRecHits.end(),RecHitLessByDet(alongMomentum));
241  outputTraj = transform(globalTrack, track, DYTRecHits);
242  }
243 
244  } else if (theMuonHitsOption == 2 ) {
245  getFirstHits(globalTrack, allRecHits, fmsRecHits);
246  outputTraj = transform(globalTrack, track, fmsRecHits);
247  }
248 
249 
250  if (outputTraj.size()) {
251  LogTrace(theCategory) << "Refitted pt: " << outputTraj.front().firstMeasurement().updatedState().globalParameters().momentum().perp() << endl;
252  return outputTraj;
253  } else {
254  LogTrace(theCategory) << "No refitted Tracks... " << endl;
255  return vector<Trajectory>();
256  }
257 
258 }
259 
260 
261 //
262 //
263 //
266  std::vector<int>& hits) const {
267 
268  LogTrace(theCategory) << " GlobalMuonRefitter::checkMuonHits " << endl;
269 
270  float coneSize = 20.0;
271  int dethits[4];
272  for ( int i=0; i<4; i++ ) hits[i]=dethits[i]=0;
273 
274  // loop through all muon hits and calculate the maximum # of hits in each chamber
275  for (ConstRecHitContainer::const_iterator imrh = all.begin(); imrh != all.end(); imrh++ ) {
276 
277  if ( (*imrh != 0 ) && !(*imrh)->isValid() ) continue;
278 
279  int station = 0;
280  int detRecHits = 0;
281  MuonRecHitContainer dRecHits;
282 
283  DetId id = (*imrh)->geographicalId();
284 
285  // Skip tracker hits
286  if (id.det()!=DetId::Muon) continue;
287 
288  if ( id.subdetId() == MuonSubdetId::DT ) {
289  DTChamberId did(id.rawId());
290  DTLayerId lid(id.rawId());
291  station = did.station();
292 
293  // Get the 1d DT RechHits from this layer
294  DTRecHitCollection::range dRecHits = theDTRecHits->get(lid);
295 
296  for (DTRecHitCollection::const_iterator ir = dRecHits.first; ir != dRecHits.second; ir++ ) {
297  double rhitDistance = fabs(ir->localPosition().x()-(**imrh).localPosition().x());
298  if ( rhitDistance < coneSize ) detRecHits++;
299  LogTrace(theCategory) << " " << (ir)->localPosition() << " " << (**imrh).localPosition()
300  << " Distance: " << rhitDistance << " recHits: " << detRecHits << endl;
301  }
302  }// end of if DT
303  else if ( id.subdetId() == MuonSubdetId::CSC ) {
304 
305  CSCDetId did(id.rawId());
306  station = did.station();
307 
308  // Get the CSC Rechits from this layer
309  CSCRecHit2DCollection::range dRecHits = theCSCRecHits->get(did);
310 
311  for (CSCRecHit2DCollection::const_iterator ir = dRecHits.first; ir != dRecHits.second; ir++ ) {
312  double rhitDistance = (ir->localPosition()-(**imrh).localPosition()).mag();
313  if ( rhitDistance < coneSize ) detRecHits++;
314  LogTrace(theCategory) << ir->localPosition() << " " << (**imrh).localPosition()
315  << " Distance: " << rhitDistance << " recHits: " << detRecHits << endl;
316  }
317  }
318  else {
319  if ( id.subdetId() != MuonSubdetId::RPC ) LogError(theCategory)<<" Wrong Hit Type ";
320  continue;
321  }
322 
323  if ( (station > 0) && (station < 5) ) {
324  if ( detRecHits > hits[station-1] ) hits[station-1] = detRecHits;
325  }
326 
327  } // end of loop over muon rechits
328 
329  for ( int i = 0; i < 4; i++ )
330  LogTrace(theCategory) <<" Station "<<i+1<<": "<<hits[i]<<" "<<dethits[i] <<endl;
331 
332  LogTrace(theCategory) << "CheckMuonHits: "<<all.size();
333 
334  // check order of muon measurements
335  if ( (all.size() > 1) &&
336  ( all.front()->globalPosition().mag() >
337  all.back()->globalPosition().mag() ) ) {
338  LogTrace(theCategory)<< "reverse order: ";
339  stable_sort(all.begin(),all.end(),RecHitLessByDet(alongMomentum));
340  }
341 }
342 
343 
344 //
345 // Get the hits from the first muon station (containing hits)
346 //
349  ConstRecHitContainer& first) const {
350 
351  LogTrace(theCategory) << " GlobalMuonRefitter::getFirstHits\nall rechits length:" << all.size() << endl;
352  first.clear();
353 
354  int station_to_keep = 999;
355  vector<int> stations;
356  for (ConstRecHitContainer::const_iterator ihit = all.begin(); ihit != all.end(); ++ihit) {
357 
358  int station = 0;
359  bool use_it = true;
360  DetId id = (*ihit)->geographicalId();
361  unsigned raw_id = id.rawId();
362  if (!(*ihit)->isValid()) station = -1;
363  else {
364  if (id.det() == DetId::Muon) {
365  switch (id.subdetId()) {
366  case MuonSubdetId::DT: station = DTChamberId(raw_id).station(); break;
367  case MuonSubdetId::CSC: station = CSCDetId(raw_id).station(); break;
368  case MuonSubdetId::RPC: station = RPCDetId(raw_id).station(); use_it = false; break;
369  }
370  }
371  }
372 
373  if (use_it && station > 0 && station < station_to_keep) station_to_keep = station;
374  stations.push_back(station);
375 
376  LogTrace(theCategory) << "rawId: " << raw_id << " station = " << station << " station_to_keep is now " << station_to_keep;
377  }
378 
379  if (station_to_keep <= 0 || station_to_keep > 4 || stations.size() != all.size())
380  LogInfo(theCategory) << "failed to getFirstHits (all muon hits are outliers/bad ?)! station_to_keep = "
381  << station_to_keep << " stations.size " << stations.size() << " all.size " << all.size();
382 
383  for (unsigned i = 0; i < stations.size(); ++i)
384  if (stations[i] >= 0 && stations[i] <= station_to_keep) first.push_back(all[i]);
385 
386  return;
387 }
388 
389 
390 //
391 // select muon hits compatible with trajectory;
392 // check hits in chambers with showers
393 //
396  const std::vector<int>& hits) const {
397 
398  ConstRecHitContainer muonRecHits;
399  const double globalChi2Cut = 200.0;
400 
401  vector<TrajectoryMeasurement> muonMeasurements = traj.measurements();
402 
403  // loop through all muon hits and skip hits with bad chi2 in chambers with high occupancy
404  for (std::vector<TrajectoryMeasurement>::const_iterator im = muonMeasurements.begin(); im != muonMeasurements.end(); im++ ) {
405 
406  if ( !(*im).recHit()->isValid() ) continue;
407  if ( (*im).recHit()->det()->geographicalId().det() != DetId::Muon ) {
408  // if ( ( chi2ndf < globalChi2Cut ) )
409  muonRecHits.push_back((*im).recHit());
410  continue;
411  }
412  ConstMuonRecHitPointer immrh = dynamic_cast<const MuonTransientTrackingRecHit*>((*im).recHit().get());
413 
414  DetId id = immrh->geographicalId();
415  int station = 0;
416  int threshold = 0;
417  double chi2Cut = 0.0;
418 
419  // get station of hit if it is in DT
420  if ( (*immrh).isDT() ) {
421  DTChamberId did(id.rawId());
422  station = did.station();
423  threshold = theHitThreshold;
424  chi2Cut = theDTChi2Cut;
425  }
426  // get station of hit if it is in CSC
427  else if ( (*immrh).isCSC() ) {
428  CSCDetId did(id.rawId());
429  station = did.station();
430  threshold = theHitThreshold;
431  chi2Cut = theCSCChi2Cut;
432  }
433  // get station of hit if it is in RPC
434  else if ( (*immrh).isRPC() ) {
435  RPCDetId rpcid(id.rawId());
436  station = rpcid.station();
437  threshold = theHitThreshold;
438  chi2Cut = theRPCChi2Cut;
439  }
440  else
441  continue;
442 
443  double chi2ndf = (*im).estimate()/(*im).recHit()->dimension();
444 
445  bool keep = true;
446  if ( (station>0) && (station<5) ) {
447  if (hits[station-1]>threshold) keep = false;
448  }
449 
450  if ( (keep || (chi2ndf<chi2Cut)) && (chi2ndf<globalChi2Cut) ) {
451  muonRecHits.push_back((*im).recHit());
452  } else {
454  << "Skip hit: " << id.det() << " " << station << ", "
455  << chi2ndf << " (" << chi2Cut << " chi2 threshold) "
456  << hits[station-1] << endl;
457  }
458  }
459 
460  // check order of rechits
461  reverse(muonRecHits.begin(),muonRecHits.end());
462  return muonRecHits;
463 }
464 
465 
466 //
467 // print RecHits
468 //
470 
471  LogTrace(theCategory) << "Used RecHits: " << hits.size();
472  for (ConstRecHitContainer::const_iterator ir = hits.begin(); ir != hits.end(); ir++ ) {
473  if ( !(*ir)->isValid() ) {
474  LogTrace(theCategory) << "invalid RecHit";
475  continue;
476  }
477 
478  const GlobalPoint& pos = (*ir)->globalPosition();
479 
481  << "r = " << sqrt(pos.x() * pos.x() + pos.y() * pos.y())
482  << " z = " << pos.z()
483  << " dimension = " << (*ir)->dimension()
484  << " det = " << (*ir)->det()->geographicalId().det()
485  << " subdet = " << (*ir)->det()->subDetector()
486  << " raw id = " << (*ir)->det()->geographicalId().rawId();
487  }
488 
489 }
490 
491 
492 //
493 // add Trajectory* to TrackCand if not already present
494 //
497 
498  if (!recHits.empty()){
499  ConstRecHitContainer::const_iterator frontHit = recHits.begin();
500  ConstRecHitContainer::const_iterator backHit = recHits.end() - 1;
501  while( !(*frontHit)->isValid() && frontHit != backHit) {frontHit++;}
502  while( !(*backHit)->isValid() && backHit != frontHit) {backHit--;}
503 
504  double rFirst = (*frontHit)->globalPosition().mag();
505  double rLast = (*backHit) ->globalPosition().mag();
506 
507  if(rFirst < rLast) return insideOut;
508  else if(rFirst > rLast) return outsideIn;
509  else {
510  LogError(theCategory) << "Impossible determine the rechits order" <<endl;
511  return undetermined;
512  }
513  } else {
514  LogError(theCategory) << "Impossible determine the rechits order" <<endl;
515  return undetermined;
516  }
517 }
518 
519 
520 //
521 // Convert Tracks into Trajectories with a given set of hits
522 //
523 vector<Trajectory> GlobalMuonRefitter::transform(const reco::Track& newTrack,
524  const reco::TransientTrack track,
525  TransientTrackingRecHit::ConstRecHitContainer recHitsForReFit) const {
526 
527  LogTrace(theCategory) << "GlobalMuonRefitter::transform: " << recHitsForReFit.size() << " hits:";
528  printHits(recHitsForReFit);
529 
530  if(recHitsForReFit.size() < 2) return vector<Trajectory>();
531 
532  // Check the order of the rechits
533  RefitDirection recHitsOrder = checkRecHitsOrdering(recHitsForReFit);
534 
535  LogTrace(theCategory) << "checkRecHitsOrdering() returned " << recHitsOrder
536  << ", theRefitDirection is " << theRefitDirection
537  << " (insideOut == " << insideOut << ", outsideIn == " << outsideIn << ")";
538 
539  // Reverse the order in the case of inconsistency between the fit direction and the rechit order
540  if(theRefitDirection != recHitsOrder) reverse(recHitsForReFit.begin(),recHitsForReFit.end());
541 
542  // Even though we checked the rechits' ordering above, we may have
543  // already flipped them elsewhere (getFirstHits() is such a
544  // culprit). Use the global positions of the states and the desired
545  // refit direction to find the starting TSOS.
546  TrajectoryStateOnSurface firstTSOS, lastTSOS;
547  unsigned int innerId, outerId;
548  bool order_swapped = track.outermostMeasurementState().globalPosition().mag() < track.innermostMeasurementState().globalPosition().mag();
549  bool inner_is_first;
550  LogTrace(theCategory) << "order swapped? " << order_swapped;
551 
552  // Fill the starting state, depending on the ordering above.
553  if ((theRefitDirection == insideOut && !order_swapped) || (theRefitDirection == outsideIn && order_swapped)) {
554  innerId = newTrack.innerDetId();
555  outerId = newTrack.outerDetId();
556  firstTSOS = track.innermostMeasurementState();
557  lastTSOS = track.outermostMeasurementState();
558  inner_is_first = true;
559  }
560  else {
561  innerId = newTrack.outerDetId();
562  outerId = newTrack.innerDetId();
563  firstTSOS = track.outermostMeasurementState();
564  lastTSOS = track.innermostMeasurementState();
565  inner_is_first = false;
566  }
567 
568  LogTrace(theCategory) << "firstTSOS: inner_is_first? " << inner_is_first
569  << " globalPosition is " << firstTSOS.globalPosition()
570  << " innerId is " << innerId;
571 
572  if(!firstTSOS.isValid()){
573  LogWarning(theCategory) << "Error wrong initial state!" << endl;
574  return vector<Trajectory>();
575  }
576 
577  firstTSOS.rescaleError(1000.);
578 
579  // This is the only way to get a TrajectorySeed with settable propagation direction
580  PTrajectoryStateOnDet garbage1;
582  PropagationDirection propDir =
583  (firstTSOS.globalPosition().basicVector().dot(firstTSOS.globalMomentum().basicVector())>0) ? alongMomentum : oppositeToMomentum;
584 
585  // These lines cause the code to ignore completely what was set
586  // above, and force propDir for tracks from collisions!
587 // if(propDir == alongMomentum && theRefitDirection == outsideIn) propDir=oppositeToMomentum;
588 // if(propDir == oppositeToMomentum && theRefitDirection == insideOut) propDir=alongMomentum;
589 
590  const TrajectoryStateOnSurface& tsosForDir = inner_is_first ? lastTSOS : firstTSOS;
591  propDir = (tsosForDir.globalPosition().basicVector().dot(tsosForDir.globalMomentum().basicVector())>0) ? alongMomentum : oppositeToMomentum;
592  LogTrace(theCategory) << "propDir based on firstTSOS x dot p is " << propDir
593  << " (alongMomentum == " << alongMomentum << ", oppositeToMomentum == " << oppositeToMomentum << ")";
594 
595  // Additional propagation diretcion determination logic for cosmic muons
596  if (theCosmicFlag) {
597  PropagationDirection propDir_first = (firstTSOS.globalPosition().basicVector().dot(firstTSOS.globalMomentum().basicVector()) > 0) ? alongMomentum : oppositeToMomentum;
598  PropagationDirection propDir_last = (lastTSOS .globalPosition().basicVector().dot(lastTSOS .globalMomentum().basicVector()) > 0) ? alongMomentum : oppositeToMomentum;
599  LogTrace(theCategory) << "propDir_first " << propDir_first << ", propdir_last " << propDir_last
600  << " : they " << (propDir_first == propDir_last ? "agree" : "disagree");
601 
602  int y_count = 0;
603  for (TransientTrackingRecHit::ConstRecHitContainer::const_iterator it = recHitsForReFit.begin(); it != recHitsForReFit.end(); ++it) {
604  if ((*it)->globalPosition().y() > 0) ++y_count;
605  else --y_count;
606  }
607 
608  PropagationDirection propDir_ycount = alongMomentum;
609  if (y_count > 0) {
610  if (theRefitDirection == insideOut) propDir_ycount = oppositeToMomentum;
611  else if (theRefitDirection == outsideIn) propDir_ycount = alongMomentum;
612  }
613  else {
614  if (theRefitDirection == insideOut) propDir_ycount = alongMomentum;
615  else if (theRefitDirection == outsideIn) propDir_ycount = oppositeToMomentum;
616  }
617 
618  LogTrace(theCategory) << "y_count = " << y_count
619  << "; based on geometrically-outermost TSOS, propDir is " << propDir << ": "
620  << (propDir == propDir_ycount ? "agrees" : "disagrees")
621  << " with ycount determination";
622 
623  if (propDir_first != propDir_last) {
624  LogTrace(theCategory) << "since first/last disagreed, using y_count propDir";
625  propDir = propDir_ycount;
626  }
627  }
628 
629  TrajectorySeed seed(garbage1,garbage2,propDir);
630 
631  if(recHitsForReFit.front()->geographicalId() != DetId(innerId)){
632  LogDebug(theCategory)<<"Propagation occured"<<endl;
633  LogTrace(theCategory) << "propagating firstTSOS at " << firstTSOS.globalPosition()
634  << " to first rechit with surface pos " << recHitsForReFit.front()->det()->surface().toGlobal(LocalPoint(0,0,0));
635  firstTSOS = theService->propagator(thePropagatorName)->propagate(firstTSOS, recHitsForReFit.front()->det()->surface());
636  if(!firstTSOS.isValid()){
637  LogDebug(theCategory)<<"Propagation error!"<<endl;
638  return vector<Trajectory>();
639  }
640  }
641 
642 /*
643  cout << " GlobalMuonRefitter : theFitter " << propDir << endl;
644  cout << " First TSOS: "
645  << firstTSOS.globalPosition() << " p="
646  << firstTSOS.globalMomentum() << " = "
647  << firstTSOS.globalMomentum().mag() << endl;
648 
649  cout << " Starting seed: "
650  << " nHits= " << seed.nHits()
651  << " tsos: "
652  << seed.startingState().parameters().position() << " p="
653  << seed.startingState().parameters().momentum() << endl;
654 
655  cout << " RecHits: "
656  << recHitsForReFit.size() << endl;
657 */
658 
659  vector<Trajectory> trajectories = theFitter->fit(seed,recHitsForReFit,firstTSOS);
660 
661  if(trajectories.empty()){
662  LogDebug(theCategory) << "No Track refitted!" << endl;
663  return vector<Trajectory>();
664  }
665 
666  return trajectories;
667 }
668 
669 
670 //
671 // Remove Selected Station Rec Hits
672 //
674 {
676  ConstRecHitContainer::const_iterator it = hits.begin();
677  for (; it!=hits.end(); it++) {
678 
679  DetId id = (*it)->geographicalId();
680 
681  //Check that this is a Muon hit that we're toying with -- else pass on this because the hacker is a moron / not careful
682 
683  if (id.det() == DetId::Tracker && theTrackerSkipSystem > 0) {
684  int layer = -999;
685  int disk = -999;
686  int wheel = -999;
687  if ( id.subdetId() == theTrackerSkipSystem){
688  // continue; //caveat that just removes the whole system from refitting
689 
690  if (theTrackerSkipSystem == PXB) {
691  PXBDetId did(id.rawId());
692  layer = did.layer();
693  }
694  if (theTrackerSkipSystem == TIB) {
695  TIBDetId did(id.rawId());
696  layer = did.layer();
697  }
698 
699  if (theTrackerSkipSystem == TOB) {
700  TOBDetId did(id.rawId());
701  layer = did.layer();
702  }
703  if (theTrackerSkipSystem == PXF) {
704  PXFDetId did(id.rawId());
705  disk = did.disk();
706  }
707  if (theTrackerSkipSystem == TID) {
708  TIDDetId did(id.rawId());
709  wheel = did.wheel();
710  }
711  if (theTrackerSkipSystem == TEC) {
712  TECDetId did(id.rawId());
713  wheel = did.wheel();
714  }
715  if (theTrackerSkipSection >= 0 && layer == theTrackerSkipSection) continue;
716  if (theTrackerSkipSection >= 0 && disk == theTrackerSkipSection) continue;
717  if (theTrackerSkipSection >= 0 && wheel == theTrackerSkipSection) continue;
718  }
719  }
720 
721  if (id.det() == DetId::Muon && theSkipStation) {
722  int station = -999;
723  int wheel = -999;
724  if ( id.subdetId() == MuonSubdetId::DT ) {
725  DTChamberId did(id.rawId());
726  station = did.station();
727  wheel = did.wheel();
728  } else if ( id.subdetId() == MuonSubdetId::CSC ) {
729  CSCDetId did(id.rawId());
730  station = did.station();
731  } else if ( id.subdetId() == MuonSubdetId::RPC ) {
732  RPCDetId rpcid(id.rawId());
733  station = rpcid.station();
734  }
735  if(station == theSkipStation) continue;
736  }
737  results.push_back(*it);
738  }
739  return results;
740 }
741 
#define LogDebug(id)
std::string theMuonRecHitBuilderName
void printHits(const ConstRecHitContainer &) const
print all RecHits of a trajectory
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
edm::Handle< DTRecHitCollection > theDTRecHits
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:53
const MuonServiceProxy * theService
unsigned int layer() const
layer id
Definition: TOBDetId.h:39
RefitDirection theRefitDirection
edm::ESHandle< TransientTrackingRecHitBuilder > theTrackerRecHitBuilder
std::string thePropagatorName
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
std::vector< Trajectory > transform(const reco::Track &newTrack, const reco::TransientTrack track, TransientTrackingRecHit::ConstRecHitContainer recHitsForReFit) const
refit the track with a new set of RecHits
T y() const
Definition: PV3DBase.h:57
GlobalPoint globalPosition() const
void setServices(const edm::EventSetup &)
set the services needed by the TrackTransformer
PropagationDirection
TransientTrackingRecHit::ConstRecHitContainer ConstRecHitContainer
ConstRecHitContainer getRidOfSelectStationHits(ConstRecHitContainer hits) const
edm::ESHandle< TransientTrackingRecHitBuilder > theMuonRecHitBuilder
unsigned int layer() const
layer id
Definition: PXBDetId.h:35
TrajectoryStateOnSurface innermostMeasurementState() const
void getFirstHits(const reco::Track &, ConstRecHitContainer &, ConstRecHitContainer &) const
get the RecHits in the tracker and the first muon chamber with hits
const edm::Event * theEvent
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:46
DataContainer const & measurements() const
Definition: Trajectory.h:169
const int keep
static const int CSC
Definition: MuonSubdetId.h:15
T mag() const
Definition: PV3DBase.h:61
edm::ESHandle< TrajectoryFitter > theFitter
T sqrt(T t)
Definition: SSEVec.h:28
double pt() const
track transverse momentum
Definition: TrackBase.h:131
T z() const
Definition: PV3DBase.h:58
MuonTransientTrackingRecHit::MuonRecHitContainer MuonRecHitContainer
unsigned int outerDetId() const
DetId of the detector on which surface the outermost state is located.
Definition: Track.h:59
edm::Handle< CSCRecHit2DCollection > theCSCRecHits
RefitDirection checkRecHitsOrdering(const ConstRecHitContainer &) const
TrajectoryStateOnSurface outermostMeasurementState() const
void checkMuonHits(const reco::Track &, ConstRecHitContainer &, std::vector< int > &) const
check muon RecHits, calculate chamber occupancy and select hits to be used in the final fit ...
TransientTrackingRecHit::ConstRecHitContainer filter(const Trajectory &)
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
bool first
Definition: L1TdeRCT.cc:79
#define LogTrace(id)
unsigned int disk() const
disk id
Definition: PXFDetId.h:43
edm::InputTag theDTRecHitLabel
std::vector< ConstRecHitPointer > ConstRecHitContainer
virtual void setEvent(const edm::Event &)
pass the Event to the algo at each event
GlobalMuonRefitter(const edm::ParameterSet &, const MuonServiceProxy *)
constructor with Parameter Set and MuonServiceProxy
Definition: DetId.h:20
const T & get() const
Definition: EventSetup.h:55
unsigned int wheel() const
wheel id
Definition: TECDetId.h:52
std::vector< Trajectory > refit(const reco::Track &globalTrack, const int theMuonHitsOption) const
build combined trajectory from sta Track and tracker RecHits
unsigned int layer() const
layer id
Definition: TIBDetId.h:41
GlobalVector globalMomentum() const
void setThr(int, int)
static const int RPC
Definition: MuonSubdetId.h:16
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
int station() const
Definition: CSCDetId.h:88
static const int DT
Definition: MuonSubdetId.h:14
virtual ~GlobalMuonRefitter()
destructor
ConstRecHitContainer selectMuonHits(const Trajectory &, const std::vector< int > &) const
select muon hits compatible with trajectory; check hits in chambers with showers
unsigned long long theCacheId_TRH
int station() const
Return the station number.
Definition: DTChamberId.h:53
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:47
T x() const
Definition: PV3DBase.h:56
std::string theTrackerRecHitBuilderName
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:54
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:61
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
edm::InputTag theCSCRecHitLabel
unsigned int wheel() const
wheel id
Definition: TIDDetId.h:50
const double par[8 *NPar][4]
T dot(const Basic3DVector &rh) const
Scalar product, or &quot;dot&quot; product, with a vector of same type.
int station() const
Definition: RPCDetId.h:98