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