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